]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSTrackerV1.cxx
Correct material budget for tracks stopped at TRD (S.Arcelli)
[u/mrichter/AliRoot.git] / ITS / AliITSTrackerV1.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15  
16 /* $Id$ */
17
18 // The purpose of this class is to permorm the ITS tracking. The 
19 // constructor has the task to inizialize some private members. The method 
20 // DoTracking is written to be called by a macro. It gets the event number,
21 // the minimum and maximum order number of TPC tracks that are to be tracked
22 // trough the ITS, and the file where the recpoints are registered. The 
23 // method Recursivetracking is a recursive function that performs the 
24 // tracking trough the ITS The method Intersection found the layer, ladder 
25 // and detector whre the intersection take place and caluclate the 
26 // cohordinates of this intersection. It returns an integer that is 0 if the 
27 // intersection has been found successfully. The two mwthods Kalmanfilter 
28 // and kalmanfiltervert operate the kalmanfilter without and with the vertex
29 // imposition respectively. The authors thank Mariana Bondila to have help 
30 // them to resolve some problems.  July-2000
31
32 #include <Riostream.h>
33 #include <Riostream.h>
34 #include <TMath.h>
35 #include <TBranch.h>
36 #include <TVector.h>
37 #include <TFile.h>
38 #include <TTree.h>
39
40 #include "TParticle.h"
41 #include "AliRun.h"
42 #include "AliITS.h"
43 #include "AliITSgeomSPD.h"
44 #include "AliITSgeomSDD.h"
45 #include "AliITSgeomSSD.h"
46 #include "AliITSgeom.h"
47 #include "AliITSRecPoint.h"
48 #include "AliKalmanTrack.h" 
49 #include "AliMagF.h"
50 #include "AliITSTrackV1.h"
51 #include "AliITSIOTrack.h"
52 #include "AliITSRad.h"   
53 #include "AliTPCtracker.h"
54 #include "AliITSTrackerV1.h"
55 #include "AliITSPid.h"
56 #include "AliMC.h"
57
58 ClassImp(AliITSTrackerV1)
59  //______________________________________________________________________
60 AliITSTrackerV1::AliITSTrackerV1() {
61   //Default constructor
62   fITS = 0;
63   fresult = 0;
64   fPtref=0.;
65   fChi2max=0.;
66   //fepsphi=0.;
67   //fepsz=0.;
68   frecPoints = 0;
69   fvettid = 0;
70   fflagvert=0;
71   frl = 0;
72   Int_t ia;
73   for(ia=0; ia<6; ia++) {
74   fNlad[ia]=0;
75   fNdet[ia]=0;
76   fAvrad[ia]=0.;
77   fDetx[ia]=0.;
78   fDetz[ia]=0.; 
79   } // end for ia  
80   fzmin = 0;
81   fzmax = 0;
82   fphimin = 0;
83   fphimax = 0;
84   fphidet = 0;
85   fNRecPoints=0;
86   fRecCylR=0;
87   fRecCylPhi=0;
88   fRecCylZ=0;
89   fFieldFactor=0;
90 }
91 //______________________________________________________________________
92 AliITSTrackerV1::AliITSTrackerV1(AliITS* IITTSS, Int_t evnumber, Bool_t flag) {
93     //Origin   A. Badala' and G.S. Pappalardo:  
94     // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
95     // Class constructor. It does some initializations.
96
97   //PH Initialisation taken from the default constructor
98     fITS      = IITTSS;
99     fresult = 0;
100     fPtref    = 0.;
101     fChi2max  =0.; 
102     frecPoints = 0;              
103     fvettid = 0;
104     fflagvert = flag;    
105     frl = 0;
106     fzmin = 0;
107     fzmax = 0;
108     fphimin = 0;
109     fphimax = 0;
110     fphidet = 0;
111   
112     Int_t imax = 200,jmax = 450;
113     frl       = new AliITSRad(imax,jmax);
114
115     //////////  gets information on geometry /////////////////////////////
116          AliITSgeom *g1 = fITS->GetITSgeom();  
117     Int_t ll=1, dd=1;
118     TVector det(9);
119
120     Int_t ia;
121     for(ia=0; ia<6; ia++) {
122         fNlad[ia]=g1->GetNladders(ia+1);
123         fNdet[ia]=g1->GetNdetectors(ia+1);
124         //cout<<fNlad[i]<<" "<<fNdet[i]<<"\n"; 
125     } // end for ia
126
127     //cout<<" mean radius = ";
128     Int_t ib;
129     for(ib=0; ib<6; ib++) {  
130         g1->GetCenterThetaPhi(ib+1,ll,dd,det);
131         Double_t r1=TMath::Sqrt(det(0)*det(0)+det(1)*det(1));
132         g1->GetCenterThetaPhi(ib+1,ll,dd+1,det);
133         Double_t r2=TMath::Sqrt(det(0)*det(0)+det(1)*det(1));
134         fAvrad[ib]=(r1+r2)/2.;
135         //cout<<fAvrad[ib]<<" ";
136     } // end for ib
137     //cout<<"\n"; getchar();
138
139     fDetx[0] = ((AliITSgeomSPD*)(g1->GetShape(1, ll, dd)))->GetDx();
140     fDetz[0] = ((AliITSgeomSPD*)(g1->GetShape(1, ll, dd)))->GetDz();
141
142     fDetx[1] = ((AliITSgeomSPD*)(g1->GetShape(2, ll, dd)))->GetDx();
143     fDetz[1] = ((AliITSgeomSPD*)(g1->GetShape(2, ll, dd)))->GetDz();
144
145     fDetx[2] = ((AliITSgeomSDD*)(g1->GetShape(3, ll, dd)))->GetDx();
146     fDetz[2] = ((AliITSgeomSDD*)(g1->GetShape(3, ll, dd)))->GetDz();
147
148     fDetx[3] = ((AliITSgeomSDD*)(g1->GetShape(4, ll, dd)))->GetDx();
149     fDetz[3] = ((AliITSgeomSDD*)(g1->GetShape(4, ll, dd)))->GetDz();
150
151     fDetx[4] = ((AliITSgeomSSD*)(g1->GetShape(5, ll, dd)))->GetDx();
152     fDetz[4] = ((AliITSgeomSSD*)(g1->GetShape(5, ll, dd)))->GetDz();
153
154     fDetx[5] = ((AliITSgeomSSD*)(g1->GetShape(6, ll, dd)))->GetDx();
155     fDetz[5] = ((AliITSgeomSSD*)(g1->GetShape(6, ll, dd)))->GetDz();
156     //cout<<"    Detx     Detz\n";
157     //for(Int_t la=0; la<6; la++) cout<<"    "<<fDetx[la]<<"     "<<
158     //                                 fDetz[la]<<endl;
159     //getchar();
160          
161     // allocate memory and define matrices fzmin, fzmax, fphimin and fphimax //
162     Double_t epsz=1.2;
163     Double_t epszdrift=0.05;
164
165     fzmin = new Double_t*[6]; fzmax = new Double_t*[6];
166     Int_t im1, im2, im2max;
167     for(im1=0; im1<6; im1++) {
168         im2max=fNdet[im1];
169         fzmin[im1] = new Double_t[im2max]; fzmax[im1] = new Double_t[im2max];
170     } // end for im1
171
172     for(im1=0; im1<6; im1++) {
173         im2max=fNdet[im1];
174         for(im2=0; im2<im2max; im2++) {
175             g1->GetCenterThetaPhi(im1+1,1,im2+1,det);
176             if(im2!=0) fzmin[im1][im2]=det(2)-fDetz[im1];
177             else   
178                 fzmin[im1][im2]=det(2)-(fDetz[im1])*epsz;
179             if(im2!=(im2max-1)) fzmax[im1][im2]=det(2)+fDetz[im1];
180             else
181                 fzmax[im1][im2]=det(2)+fDetz[im1]*epsz;
182             if(im1==2 || im1==3) {
183                 fzmin[im1][im2]-=epszdrift;
184                 fzmax[im1][im2]+=epszdrift;
185             } // end if im1==2 || im1==3
186         } // end for im2
187     } // end for im1
188
189     fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
190     for(im1=0;im1<6;im1++) {
191         im2max=fNlad[im1];
192         fphimin[im1] = new Double_t[im2max]; 
193         fphimax[im1] = new Double_t[im2max];
194     } // end for im1
195   
196     fphidet = new Double_t*[6];
197     for(im1=0; im1<6; im1++) {
198         im2max=fNlad[im1];
199         fphidet[im1] = new Double_t[im2max];
200     } // end for im1
201
202         Double_t global[3],local[3];
203     Double_t pigre=TMath::Pi();
204     Double_t xmin,ymin,xmax,ymax;
205
206     for(im1=0; im1<6; im1++) {
207         im2max=fNlad[im1];
208         for(im2=0; im2<im2max; im2++) {
209             Int_t idet=2;
210             g1->GetCenterThetaPhi(im1+1,im2+1,idet,det);
211             fphidet[im1][im2] = TMath::ATan2(Double_t(det(1)),
212                                              Double_t(det(0))); 
213             if(fphidet[im1][im2]<0.) fphidet[im1][im2]+=2.*pigre;  
214             local[1]=local[2]=0.;  
215             local[0]= -(fDetx[im1]);    
216             if(im1==0) local[0]= (fDetx[im1]); //to take into account 
217                                                // different reference system
218             g1->LtoG(im1+1,im2+1,idet,local,global);
219             xmax=global[0]; ymax=global[1];
220             local[0]= (fDetx[im1]);   
221             if(im1==0) local[0]= -(fDetx[im1]);//take into account different 
222                                                // reference system
223             g1->LtoG(im1+1,im2+1,idet,local,global);
224             xmin=global[0]; ymin=global[1];
225             fphimin[im1][im2]= TMath::ATan2(ymin,xmin);
226             if(fphimin[im1][im2]<0.) fphimin[im1][im2]+=2.*pigre; 
227             fphimax[im1][im2]= TMath::ATan2(ymax,xmax);
228             if(fphimax[im1][im2]<0.) fphimax[im1][im2]+=2.*pigre;
229         } // end for im2
230     } // end for im1
231 //////////////////////////////////////////////////////////////////////////////////////////////////////////
232 /////////////// allocate memory and define vector fNRecPoints and matrices fRecCylR, fRecCylPhi, fRecCylZ /////////////
233         gAlice->GetEvent(evnumber);
234   Int_t numOfModules = g1->GetIndexMax();
235   fRecCylR = new Double_t *[numOfModules];
236   fRecCylPhi = new Double_t *[numOfModules]; 
237   fRecCylZ = new Double_t *[numOfModules];
238   AliITSRecPoint *recp;
239   fNRecPoints = new Int_t[numOfModules];
240    
241                  for(Int_t module=0; module<numOfModules; module++) {                           
242                   fITS->ResetRecPoints();                    
243         gAlice->TreeR()->GetEvent(module);                
244                   frecPoints=fITS->RecPoints();
245                   Int_t nRecPoints=fNRecPoints[module]=frecPoints->GetEntries();
246                   fRecCylR[module] = new Double_t[nRecPoints];
247                   fRecCylPhi[module] = new Double_t[nRecPoints];
248                   fRecCylZ[module] = new  Double_t[nRecPoints];           
249                   Int_t ind;
250                   for(ind=0; ind<fNRecPoints[module]; ind++) {    
251                     recp=(AliITSRecPoint*)frecPoints->UncheckedAt(ind);                                          
252                         // Float_t global[3], local[3];
253                          Double_t global[3], local[3];
254                local[0]=recp->GetX();
255                local[1]=0.;
256                local[2]= recp->GetZ();                                  
257                          g1->LtoG(module,local,global);
258                                                                                          
259                          Double_t r = TMath::Sqrt(global[0]*global[0]+global[1]*global[1]);                     // r hit
260                          Double_t phi = TMath::ATan2(global[1],global[0]); if(phi<0.) phi+=2.*TMath::Pi();      // phi hit                      
261                  Double_t z = global[2];                                                                // z hit
262                                                                                                                                                                                                                                      
263                          fRecCylR[module][ind]=r;
264                          fRecCylPhi[module][ind]=phi;
265                          fRecCylZ[module][ind]=z;                        
266                   }             
267                 }        
268          //}  
269   //}
270 ////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
271          
272
273     ////////// gets magnetic field factor //////////////////////////////
274
275     AliMagF * fieldPointer = gAlice->Field();
276    // fFieldFactor = (Double_t)fieldPointer->Factor();
277     fFieldFactor =(Double_t)fieldPointer-> SolenoidField()/10/.2;
278    // cout<< " field factor = "<<fFieldFactor<<"\n"; getchar();
279 }
280 //______________________________________________________________________
281 AliITSTrackerV1::AliITSTrackerV1(const AliITSTrackerV1 &cobj) : TObject(cobj) {
282     // Origin  A. Badala' and G.S. Pappalardo:
283     // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
284     // copy constructor
285          
286          *fITS = *cobj.fITS;
287     *fresult = *cobj.fresult;
288     fPtref = cobj.fPtref;
289     fChi2max = cobj.fChi2max;    
290     **fvettid = **cobj.fvettid;
291     fflagvert = cobj.fflagvert;
292     Int_t imax=200,jmax=450;
293     frl = new AliITSRad(imax,jmax);      
294     *frl = *cobj.frl;
295     fFieldFactor = cobj.fFieldFactor;
296     Int_t i,im1,im2,im2max;
297     for(i=0; i<6; i++) {
298         fNlad[i] = cobj.fNlad[i];
299         fNdet[i] = cobj.fNdet[i]; 
300         fAvrad[i] = cobj.fAvrad[i];
301         fDetx[i] = cobj.fDetx[i];
302         fDetz[i] = cobj.fDetz[i];
303     } // end or i
304     fzmin = new Double_t*[6]; fzmax = new Double_t*[6];
305     for(im1=0; im1<6; im1++) {
306         im2max=fNdet[im1];
307         fzmin[im1] = new Double_t[im2max];
308         fzmax[im1] = new Double_t[im2max];
309     } // end for im1
310     fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
311     for(im1=0;im1<6;im1++) {
312         im2max=fNlad[im1];
313         fphimin[im1] = new Double_t[im2max];
314         fphimax[im1] = new Double_t[im2max];
315     } // end for im1
316   
317     fphidet = new Double_t*[6];
318     for(im1=0; im1<6; im1++) {
319         im2max=fNlad[im1];
320         fphidet[im1] = new Double_t[im2max];
321     } // end for im1
322     for(im1=0; im1<6; im1++) {
323         im2max=fNdet[im1];
324         for(im2=0; im2<im2max; im2++) {
325             fzmin[im1][im2]=cobj.fzmin[im1][im2];
326             fzmax[im1][im2]=cobj.fzmax[im1][im2];
327         } // end for im2
328     } // end for im1
329     for(im1=0; im1<6; im1++) {
330         im2max=fNlad[im1];
331         for(im2=0; im2<im2max; im2++) {
332             fphimin[im1][im2]=cobj.fphimin[im1][im2];
333             fphimax[im1][im2]=cobj.fphimax[im1][im2];
334             fphidet[im1][im2]=cobj.fphidet[im1][im2];  
335         } // end for im2
336     } // end for im2
337
338
339         AliITSgeom *g1 = fITS->GetITSgeom();  
340    Int_t numOfModules = g1->GetIndexMax();
341         /*
342   fRecCylR = new Float_t *[numOfModules];
343   fRecCylPhi = new Float_t *[numOfModules]; 
344   fRecCylZ = new Float_t *[numOfModules];
345   */
346   fRecCylR = new Double_t *[numOfModules];
347   fRecCylPhi = new Double_t *[numOfModules]; 
348   fRecCylZ = new Double_t *[numOfModules];  
349   fNRecPoints = new Int_t[numOfModules];        
350                 for(Int_t module=0; module<numOfModules; module++) {            
351                   Int_t nRecPoints=fNRecPoints[module]=cobj.fNRecPoints[module];
352                   /*
353                   fRecCylR[module] = new Float_t[nRecPoints];
354                   fRecCylPhi[module] = new Float_t[nRecPoints];
355                   fRecCylZ[module] = new Float_t[nRecPoints];
356                   */
357                   fRecCylR[module] = new Double_t[nRecPoints];
358                   fRecCylPhi[module] = new Double_t[nRecPoints];
359                   fRecCylZ[module] = new Double_t[nRecPoints];            
360                   Int_t ind;    
361                   for(ind=0; ind<nRecPoints; ind++) {       
362                          fRecCylR[module][ind]=cobj.fRecCylR[module][ind];
363                          fRecCylPhi[module][ind]=cobj.fRecCylPhi[module][ind];
364                          fRecCylZ[module][ind]=cobj.fRecCylZ[module][ind];                       
365                   }             
366                 }        
367  
368 }
369 /*
370 //______________________________________________________________________
371 void AliITSTrackerV1::DelMatrix(Int_t numOfModules) {
372   // cleanup of some data members
373   for(Int_t mod=0; mod<numOfModules; mod++) {
374     delete fRecCylR[mod];
375          delete fRecCylPhi[mod];
376          delete fRecCylZ[mod];
377   }
378     delete fRecCylR;
379          delete fRecCylPhi;
380          delete fRecCylZ;
381 }
382 */
383 //______________________________________________________________________
384 AliITSTrackerV1::~AliITSTrackerV1(){
385     // Origin  A. Badala' and G.S. Pappalardo:
386     // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it  
387     // class destructor  
388   delete frl;
389   delete fNRecPoints;
390   for(Int_t i=0; i<6; i++) {
391     delete fzmin[i];
392          delete fzmax[i];
393          delete fphimin[i];
394          delete fphimax[i];
395          delete fphidet[i];
396   }
397
398   delete fzmin;
399   delete fzmax;
400   delete fphimin;
401   delete fphimax;
402   delete fphidet;
403          
404 }
405 //______________________________________________________________________
406 AliITSTrackerV1 &AliITSTrackerV1::operator=(const AliITSTrackerV1 &obj) {
407     // Origin  A. Badala' and G.S. Pappalardo:
408     // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it  
409     // assignement operator
410
411          *fITS = *obj.fITS;
412     *fresult = *obj.fresult;
413     fPtref = obj.fPtref;
414     fChi2max = obj.fChi2max;      
415     **fvettid = **obj.fvettid;
416     fflagvert = obj.fflagvert;
417     Int_t imax=200,jmax=450;
418     frl = new AliITSRad(imax,jmax);      
419     *frl = *obj.frl;
420     fFieldFactor = obj.fFieldFactor;
421     Int_t i;
422     for(i=0; i<6; i++) {
423         fNlad[i] = obj.fNlad[i];
424         fNdet[i] = obj.fNdet[i]; 
425         fAvrad[i] = obj.fAvrad[i];
426         fDetx[i] = obj.fDetx[i];
427         fDetz[i] = obj.fDetz[i];
428     } // end for i
429     fzmin = new Double_t*[6];
430     fzmax = new Double_t*[6];
431     Int_t im1, im2, im2max;
432     for(im1=0; im1<6; im1++) {
433         im2max=fNdet[im1];
434         fzmin[im1] = new Double_t[im2max]; fzmax[im1] = new Double_t[im2max];
435     } // end for im1
436     fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
437     for(im1=0;im1<6;im1++) {
438         im2max=fNlad[im1];
439         fphimin[im1] = new Double_t[im2max];
440         fphimax[im1] = new Double_t[im2max];
441     } // end for im1
442
443     fphidet = new Double_t*[6];
444     for(im1=0; im1<6; im1++) {
445         im2max=fNlad[im1];
446         fphidet[im1] = new Double_t[im2max];
447     } // end for im1
448     for(im1=0; im1<6; im1++) {
449         im2max=fNdet[im1];
450         for(im2=0; im2<im2max; im2++) {
451             fzmin[im1][im2]=obj.fzmin[im1][im2];
452             fzmax[im1][im2]=obj.fzmax[im1][im2];
453         } // end for im2
454     } // end for im1
455     for(im1=0; im1<6; im1++) {
456         im2max=fNlad[im1];
457         for(im2=0; im2<im2max; im2++) {
458             fphimin[im1][im2]=obj.fphimin[im1][im2];
459             fphimax[im1][im2]=obj.fphimax[im1][im2];
460             fphidet[im1][im2]=obj.fphidet[im1][im2];  
461         } // end for im2
462     } // end for im1
463
464         AliITSgeom *g1 = fITS->GetITSgeom();  
465    Int_t numOfModules = g1->GetIndexMax();
466   fRecCylR = new Double_t *[numOfModules];
467   fRecCylPhi = new Double_t *[numOfModules]; 
468   fRecCylZ = new Double_t *[numOfModules];  
469   fNRecPoints = new Int_t[numOfModules];  
470           for(Int_t module=0; module<numOfModules; module++) {            
471                   Int_t nRecPoints=fNRecPoints[module]=obj.fNRecPoints[module];
472                   fRecCylR[module] = new Double_t[nRecPoints];
473                   fRecCylPhi[module] = new Double_t[nRecPoints];
474                   fRecCylZ[module] = new Double_t[nRecPoints];            
475                   Int_t ind;
476                   for(ind=0; ind<nRecPoints; ind++) {     
477                          fRecCylR[module][ind]=obj.fRecCylR[module][ind];
478                          fRecCylPhi[module][ind]=obj.fRecCylPhi[module][ind];
479                          fRecCylZ[module][ind]=obj.fRecCylZ[module][ind];                        
480                   }             
481                 }        
482          
483          
484     return *this;
485 }
486 //______________________________________________________________________
487 void AliITSTrackerV1::DoTracking(Int_t evNumber,Int_t minTr,Int_t maxTr,
488                                  TFile *file, Bool_t realmass) {
489     // Origin   A. Badala' and G.S. Pappalardo:
490     // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
491     // The method needs the event number, the minimum and maximum order
492     // number of TPC tracks that 
493     // are to be tracked trough the ITS, and the file where the recpoints
494     // are registered.
495     // The method can be called by a macro. It preforms the tracking for
496     // all good TPC tracks
497
498     printf("begin DoTracking - file %p\n",(void*)file);
499
500     gAlice->GetEvent(evNumber);  //modificato per gestire hbt
501  
502     AliKalmanTrack::SetConvConst(1000/0.299792458/gAlice->Field()->SolenoidField());
503    // cout<<" field = "<<gAlice->Field()->SolenoidField()<<endl;
504
505
506     TFile *cf=TFile::Open("AliTPCclusters.root");  
507     AliTPCParam *digp= (AliTPCParam*)cf->Get("75x40_100x60_150x60");
508     if (!digp) { cerr<<"TPC parameters have not been found !\n"; getchar();}
509
510     cf->cd();
511     TString foldname(fITS->GetLoader()->GetEventFolder()->GetName());
512     
513     printf("This method is not converted to the NewIO !\n");  //I.B.
514     return; //I.B.
515     AliTPCtracker *tracker = new AliTPCtracker(digp);            //I.B.
516     //PH    tracker->SetEventNumber(evNumber);                           //I.B.
517
518     // Load clusters
519     printf("This method is not converted to the NewIO !\n");  //I.B.
520     return; //I.B.
521     tracker->LoadClusters(0);       //I.B.
522
523     // Load tracks
524     TFile *tf=TFile::Open("AliTPCtracksSorted.root");  
525     if (!tf->IsOpen()) {
526         cerr<<"Can't open AliTPCtracksSorted.root !\n";
527         return ;
528     } // end if
529     TObjArray tracks(200000);
530     char tname[100];                 
531     sprintf(tname,"TreeT_TPC_%d",evNumber);
532
533     TTree *tracktree=(TTree*)tf->Get(tname);   
534     if (!tracktree) {cerr<<"Can't get a tree with TPC tracks !\n";}   
535     TBranch *tbranch=tracktree->GetBranch("tracks");
536     Int_t nentr=(Int_t)tracktree->GetEntries();
537     Int_t kk;
538
539     AliITSRecPoint *recp;    // oggi
540         AliTPCtrack *ioTrackTPC=0;    
541     for (kk=0; kk<nentr; kk++) {
542         ioTrackTPC=new AliTPCtrack; 
543         tbranch->SetAddress(&ioTrackTPC);
544         tracktree->GetEvent(kk);    
545         tracker->CookLabel(ioTrackTPC,0.1);       
546         tracks.AddLast(ioTrackTPC);         
547     } // end for kk
548     delete tracker;      
549     tf->Close();
550
551     Int_t nt = tracks.GetEntriesFast();
552     cerr<<"Number of found tracks "<<nt<<endl;
553
554     TVector vec(5);
555     TTree *tr=gAlice->TreeR();
556     Int_t nent=(Int_t)tr->GetEntries();  
557     frecPoints = fITS->RecPoints();
558
559     Int_t numbpoints;
560     Int_t totalpoints=0;
561     Int_t *np = new Int_t[nent];
562     fvettid = new Int_t* [nent];
563     Int_t mod;
564
565     for (mod=0; mod<nent; mod++) {
566         fvettid[mod]=0;
567         fITS->ResetRecPoints();  
568         gAlice->TreeR()->GetEvent(mod); 
569         numbpoints = frecPoints->GetEntries();
570         totalpoints+=numbpoints;
571         np[mod] = numbpoints;
572         //cout<<" mod = "<<mod<<"   numbpoints = "<<numbpoints<<"\n";getchar();
573         fvettid[mod] = new Int_t[numbpoints];
574         Int_t ii;
575         for (ii=0;ii<numbpoints; ii++) *(fvettid[mod]+ii)=0;
576     } // end for mod
577
578     AliTPCtrack *track=0;
579
580     if(minTr < 0) {minTr = 0; maxTr = nt-1;}   
581
582     TVector vgeant(3);
583
584     TTree tracktree1("TreeT","Tree with ITS tracks");
585     AliITSIOTrack *ioTrack=0;
586         AliITSPid *pid=new AliITSPid(1000);  // oggi
587
588     tracktree1.Branch("ITStracks","AliITSIOTrack",&ioTrack,32000,0);
589   
590     TDatabasePDG * db = new TDatabasePDG;   
591   
592     Int_t j;       
593     for (j=minTr; j<=maxTr; j++) {     
594         track=(AliTPCtrack*)tracks.UncheckedAt(j);
595         if (!track) continue;
596         
597         ///   mass definition ////////////////////////
598     Double_t mass=0.13956995;
599          Int_t pcode=211;  // a pion by default
600          
601          if(realmass) {
602         if(TMath::Abs(pcode)<20443) mass=db->GetParticle(pcode)->Mass();
603          }
604          else {
605                  mass = track->GetMass();
606 //               cout << "Mass = " << mass << endl;
607          }
608          
609         
610         
611          //   new propagation to the end of TPC
612     Double_t xk=80.;
613    // track->PropagateTo(xk,0.,0.); //Ne if it's still there   //attenzione funziona solo se modifica in TPC
614         // Double_t xk=77.415;   
615          track->PropagateTo(xk, 28.94, 1.204e-3);
616     xk-=0.005;
617     track->PropagateTo(xk, 44.77,1.71); //Tedlar         
618     xk-=0.02;
619     track->PropagateTo(xk, 44.86, 1.45);   //Kevlar
620     xk-=2.0;
621     track->PropagateTo(xk, 41.28, 0.029);//Nomex
622     xk-=0.02;
623     track->PropagateTo(xk, 44.86, 1.45);   //Kevlar
624     xk-=0.005;
625     track->PropagateTo(xk, 44.77, 1.71); //Tedlar
626
627     xk=61.;
628    // track->PropagateTo(xk,0.,0.); //C02
629          track->PropagateTo(xk,36.2,1.98e-3); //C02        //attenzione funziona solo se modifica in TPC
630
631     xk -=0.005;
632     track->PropagateTo(xk, 24.01, 2.7);    //Al    
633     xk -=0.005;
634     track->PropagateTo(xk, 44.77, 1.71);  //Tedlar
635     xk -=0.02;
636     track->PropagateTo(xk, 44.86, 1.45);    //Kevlar
637     xk -=0.5;
638     track->PropagateTo(xk, 41.28, 0.029);  //Nomex    
639     xk -=0.02;
640     track->PropagateTo(xk, 44.86, 1.45);    //Kevlar
641     xk -=0.005;
642     track->PropagateTo(xk, 44.77, 1.71);  //Tedlar
643     xk -=0.005;
644     track->PropagateTo(xk, 24.01, 2.7);    //Al 
645     
646         ////////////////////////////////////////////////////////////////////////////////////////////////////////        
647         //AliITSTrackV1 trackITS(*track);
648         AliITSTrackV1 trackITS(*track, fFieldFactor);
649         //cout<<" fFieldFactor = "<<fFieldFactor<<"\n";      
650         trackITS.PutMass(mass);   //new to add mass to track
651         if(fresult){ delete fresult; fresult=0;}         
652         fresult = new AliITSTrackV1(trackITS);   
653
654         AliITSTrackV1 primaryTrack(trackITS);
655         vgeant=(*fresult).GetVertex();
656
657         // Definition of dv and zv for vertex constraint        
658         Double_t sigmaDv=0.0050;  Double_t sigmaZv=0.010;       
659         //Double_t sigmaDv=0.0015;  Double_t sigmaZv=0.0015;
660         Double_t uniform= gRandom->Uniform();
661         Double_t signdv;
662         if(uniform<=0.5) signdv=-1.;
663         else
664             signdv=1.;
665
666         Double_t vr=TMath::Sqrt(vgeant(0)*vgeant(0)+ vgeant(1)*vgeant(1));
667         Double_t dv=gRandom->Gaus(signdv*vr,(Float_t)sigmaDv); 
668         Double_t zv=gRandom->Gaus(vgeant(2),(Float_t)sigmaZv);          
669         //cout<<" Dv e Zv = "<<dv<<" "<<zv<<"\n";
670         trackITS.SetDv(dv);
671         trackITS.SetZv(zv);
672         trackITS.SetsigmaDv(sigmaDv);
673         trackITS.SetsigmaZv(sigmaZv); 
674         (*fresult).SetDv(dv);
675         (*fresult).SetZv(zv);
676         (*fresult).SetsigmaDv(sigmaDv);
677         (*fresult).SetsigmaZv(sigmaZv);
678         primaryTrack.SetDv(dv);
679         primaryTrack.SetZv(zv);
680         primaryTrack.SetsigmaDv(sigmaDv);
681         primaryTrack.SetsigmaZv(sigmaZv);
682         primaryTrack.PrimaryTrack();
683         TVector  d2=primaryTrack.Getd2();
684         TVector  tgl2=primaryTrack.Gettgl2();
685         TVector  dtgl=primaryTrack.Getdtgl();
686         trackITS.Setd2(d2); trackITS.Settgl2(tgl2);
687         trackITS.Setdtgl(dtgl); 
688         (*fresult).Setd2(d2); (*fresult).Settgl2(tgl2);
689         (*fresult).Setdtgl(dtgl);
690         /*
691           trackITS.SetVertex(vertex); trackITS.SetErrorVertex(ervertex);
692           (*result).SetVertex(vertex);   (*result).SetErrorVertex(ervertex);   
693         */
694         TList *list= new TList();
695
696         list->AddLast(&trackITS);
697   
698         fPtref=TMath::Abs( (trackITS).GetPt() );
699         //cout<<" fPtref = " <<fPtref<<"\n";
700         if(fPtref>1.0) fChi2max=40.;         
701         if(fPtref<=1.0) fChi2max=20.;
702         if(fPtref<0.4 ) fChi2max=100.;
703         if(fPtref<0.2 ) fChi2max=40.;             
704         // if(fPtref<0.4 ) fChi2max=30.;                                 
705         // if(fPtref<0.2 ) fChi2max=20.;
706         //if(fPtref<0.2 ) fChi2max=10.;
707         //if(fPtref<0.1 ) fChi2max=5.;
708         //cout << "\n Pt = " << fPtref <<"\n";  //stampa
709         RecursiveTracking(list);   
710         list->Delete();
711         delete list;
712
713         Int_t itot=-1;
714         TVector vecTotLabRef(18);
715         Int_t lay, k;
716         for(lay=5; lay>=0; lay--) {
717             TVector vecLabRef(3); 
718             vecLabRef=(*fresult).GetLabTrack(lay);
719             Float_t clustZ=(*fresult).GetZclusterTrack( lay);   
720             for(k=0; k<3; k++){  
721                 Int_t lpp=(Int_t)vecLabRef(k);
722                 if(lpp>=0) {
723                     TParticle *p=(TParticle*) gAlice->GetMCApp()->Particle(lpp);
724                     Int_t pcode=p->GetPdgCode();
725                     if(pcode==11) vecLabRef(k)=p->GetFirstMother();
726                 } // end if
727                 itot++; vecTotLabRef(itot)=vecLabRef(k);
728                 if(vecLabRef(k)==0. && clustZ == -1.) vecTotLabRef(itot) =-3.;
729             } // end for k
730         } // end for lay
731         Long_t labref;
732         Int_t freq;  
733         (*fresult).Search(vecTotLabRef, labref, freq);
734
735         //if(freq < 6) labref=-labref;        // cinque - sei
736         if(freq < 5) labref=-labref;        // cinque - sei     
737         (*fresult).SetLabel(labref);
738
739         // cout<<" progressive track number = "<<j<<"\r";
740         // cout<<j<<"\r";
741         Int_t numOfCluster=(*fresult).GetNumClust();  
742         //cout<<" progressive track number = "<<j<<"\n";    // stampa
743         Long_t labITS=(*fresult).GetLabel();
744         //cout << " ITS track label = " << labITS << "\n";      // stampa           
745         Int_t lab=track->GetLabel();                
746         //cout << " TPC track label = " << lab <<"\n";      // stampa
747         //propagation to vertex
748
749         Double_t rbeam=3.;
750         if((*fresult).DoNotCross(rbeam)) continue;  //no intersection with beampipe     
751         (*fresult).Propagation(rbeam);
752         Double_t c00,c10,c11,c20,c21,c22,c30,c31,c32,c33,c40,c41,c42,c43,c44;
753         (*fresult).GetCElements(c00,
754                                 c10,c11,
755                                 c20,c21,c22,
756                                 c30,c31,c32,c33,
757                                 c40,c41,c42,c43,c44);
758                  
759         Double_t pt=TMath::Abs((*fresult).GetPt());
760         Double_t dr=(*fresult).GetD();
761         Double_t z=(*fresult).GetZ();
762         Double_t tgl=(*fresult).GetTgl();
763         Double_t c=(*fresult).GetC();
764         Double_t cy=c/2.;
765         Double_t dz=z-(tgl/cy)*TMath::ASin((*fresult).Arga(rbeam));
766         dz-=vgeant(2);
767         // cout<<" dr e dz alla fine = "<<dr<<" "<<dz<<"\n"; getchar();
768         Double_t phi=(*fresult).Getphi();
769         Double_t phivertex = phi - TMath::ASin((*fresult).ArgA(rbeam));
770         Double_t duepi=2.*TMath::Pi();   
771         if(phivertex>duepi) phivertex-=duepi;
772         if(phivertex<0.) phivertex+=duepi;
773         /////////////////////////////////////////////////////////////
774         Int_t idmodule,idpoint;
775         if(numOfCluster >=5)  {            // cinque - sei
776             //if(numOfCluster ==6)  {            // cinque - sei
777             AliITSIOTrack outTrack;
778             ioTrack=&outTrack;
779             ioTrack->SetStatePhi(phi);
780             ioTrack->SetStateZ(z);
781             ioTrack->SetStateD(dr);
782             ioTrack->SetStateTgl(tgl);
783             ioTrack->SetStateC(c);
784             Double_t radius=(*fresult).Getrtrack();
785             ioTrack->SetRadius(radius);
786             Int_t charge;
787             if(c>0.) charge=-1;  else charge=1;
788             ioTrack->SetCharge(charge);
789                 Double_t trackmass=(*fresult).GetMass();         // oggi
790                 ioTrack->SetMass(trackmass);                     // oggi
791             ioTrack->SetCovMatrix(c00,
792                                   c10,c11,
793                                   c20,c21,c22,
794                                   c30,c31,c32,c33,
795                                   c40,c41,c42,c43,c44);
796             Double_t px=pt*TMath::Cos(phivertex);
797             Double_t py=pt*TMath::Sin(phivertex);
798             Double_t pz=pt*tgl;
799             Double_t xtrack=dr*TMath::Sin(phivertex);
800             Double_t ytrack=dr*TMath::Cos(phivertex);
801             Double_t ztrack=dz+vgeant(2);
802             ioTrack->SetPx(px);
803             ioTrack->SetPy(py);
804             ioTrack->SetPz(pz);
805             ioTrack->SetX(xtrack);
806             ioTrack->SetY(ytrack);
807             ioTrack->SetZ(ztrack);
808             ioTrack->SetLabel(labITS);
809             ioTrack->SetTPCLabel(lab);
810                  ioTrack->SetDz(dz);
811             Int_t il;           
812             /*
813                 for(il=0;il<6; il++){
814                 ioTrack->SetIdPoint(il,(*fresult).GetIdPoint(il));
815                 ioTrack->SetIdModule(il,(*fresult).GetIdModule(il));
816             } // end for il
817             */
818             //tracktree1.Fill();
819                 Float_t q[4]={-1.,-1.,-1.,-1.};
820             Float_t  globaldedx=0.;        
821             for (il=0;il<6;il++) {
822                 idpoint=(*fresult).GetIdPoint(il);
823                 idmodule=(*fresult).GetIdModule(il);
824                 if(idmodule>0.) *(fvettid[idmodule]+idpoint)=1;
825                                                                
826                 ioTrack->SetIdPoint(il,idpoint);
827                 ioTrack->SetIdModule(il,idmodule);
828                 ////  for q definition
829                 if(il>1){
830                   if(idmodule>0.){                      
831                     fITS->ResetRecPoints();
832                     gAlice->TreeR()->GetEvent(idmodule);
833                     recp=(AliITSRecPoint*)frecPoints->UncheckedAt(idpoint);
834                     q[il-2]=recp->GetQ()*(*fresult).Getfcor(il-2);
835                   }
836                 }                               
837             } // end for il      
838                 q[0]/=280.; q[1]/=280.;
839             q[2]/=38.; q[3]/=38.;
840
841  // cout<<" q prima = "<<q[0]<<" "<<q[1]<<" "<<q[2]<<" "<<q[3]<<"\n"; getchar(); 
842      
843      Int_t swap;
844   do{
845     swap=0;   
846     for (il=0; il<3; il++) {
847       if (q[il]<=q[il+1]) continue;
848       Float_t tmp=q[il];
849       q[il]=q[il+1]; q[il+1]=tmp;
850       swap++;
851     }
852   } while(swap); 
853
854  
855   // cout<<" q dopo = "<<q[0]<<" "<<q[1]<<" "<<q[2]<<" "<<q[3]<<"\n"; getchar();
856       
857     if(q[0]<0.) {
858       q[0]=q[1];
859       q[1]=q[2];
860       q[2]=q[3];
861       q[3]=-1.;             
862     } 
863  
864   // cout<<" q  dopo if = "<<q[0]<<" "<<q[1]<<" "<<q[2]<<" "<<q[3]<<"\n"; getchar(); 
865      
866     globaldedx=(q[0]+q[1])/2.;
867     
868    // if(q[3]> 0.) globaldedx=(q[0]+q[1]+q[2]+q[3])/4.;
869    //      else    globaldedx=(q[0]+q[1]+q[2])/3.; 
870    
871     ioTrack->SetdEdx(globaldedx);
872     ioTrack->SetPid(pid->GetPcode(ioTrack));
873             
874             tracktree1.Fill();         
875         } // end if on numOfCluster
876         //gObjectTable->Print();    // stampa memoria     
877     }  //  end for (int j=minTr; j<=maxTr; j++)
878     delete db;           
879     static Bool_t first=kTRUE;
880     static TFile *tfile;
881     if(first) {
882         tfile=new TFile("itstracks.root","RECREATE");
883         //cout<<"I have opened itstracks.root file "<<endl;
884     } // end if     
885     first=kFALSE;
886     tfile->cd();
887     tfile->ls();
888     char hname[30];
889     sprintf(hname,"TreeT%d",evNumber);
890         cout << "Number of saved ITS tracks " << tracktree1.GetEntries() << endl;
891     tracktree1.Write(hname);
892   
893     TTree *fAli=gAlice->TreeK();
894     TFile *fileAli=0;
895     if (fAli) fileAli =fAli->GetCurrentFile();
896     fileAli->cd();
897     ////////////////////////////////////////////////////////////////////
898
899     printf("delete vectors\n");
900     if(np) delete [] np;
901     if(fvettid) delete [] fvettid;
902     if(fresult) {delete fresult; fresult=0;}
903 }
904 //______________________________________________________________________
905 void AliITSTrackerV1::RecursiveTracking(TList *trackITSlist) {
906     //   This function perform the recursive tracking in ITS detectors
907     // reference is a pointer to the final best track
908     // Origin  A. Badala' and G.S. Pappalardo:
909     // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
910     // The authors thank Mariana Bondila to have help them to resolve some
911     // problems.  July-2000
912
913     //Rlayer[0]=4.; Rlayer[1]=7.;  Rlayer[2]=14.9;
914     // Rlayer[3]=23.8;  Rlayer[4]=39.1;  Rlayer[5]=43.6; //vecchio
915
916     ////////////////////// 
917     Float_t sigmaphil[6], sigmazl[6];
918     sigmaphil[0]=1.44e-6/(fAvrad[0]*fAvrad[0]);
919     sigmaphil[1]=1.44e-6/(fAvrad[1]*fAvrad[1]);
920     sigmaphil[2]=1.444e-5/(fAvrad[2]*fAvrad[2]);
921     sigmaphil[3]=1.444e-5/(fAvrad[3]*fAvrad[3]);
922     sigmaphil[4]=4e-6/(fAvrad[4]*fAvrad[4]);
923     sigmaphil[5]=4e-6/(fAvrad[5]*fAvrad[5]);
924     sigmazl[0]=1e-2;
925     sigmazl[1]=1e-2;
926     sigmazl[2]=7.84e-4;
927     sigmazl[3]=7.84e-4;
928     sigmazl[4]=0.6889;
929     sigmazl[5]=0.6889;  
930     ///////////////////////////////////////////////////////////
931     Int_t index; 
932     AliITSgeom *g1 = fITS->GetITSgeom();
933     AliITSRecPoint *recp;        
934     for(index =0; index<trackITSlist->GetSize(); index++) {
935         AliITSTrackV1 *trackITS = (AliITSTrackV1 *) trackITSlist->At(index);
936         if((*trackITS).GetLayer()==7) fresult->SetChi2(10.223e140);
937         // cout <<" Layer inizio = "<<(*trackITS).GetLayer()<<"\n";
938         //  cout<<"fvtrack =" <<"\n";
939         //  cout << (*trackITS)(0) << " "<<(*trackITS)(1)<<" "
940         //       <<(*trackITS)(2)<<" "<<(*trackITS)(3)<<" "
941         //       <<(*trackITS)(4)<<"\n";
942         //  cout<< " rtrack = "<<(*trackITS).Getrtrack()<<"\n";
943         //  cout<< " Pt = "<<(*trackITS).GetPt()<<"\n";
944         //  getchar();    
945         Double_t chi2Now, chi2Ref;
946         Float_t numClustRef = fresult->GetNumClust();            
947         if((*trackITS).GetLayer()==1 ) {
948             chi2Now = trackITS->GetChi2();
949             Float_t numClustNow = trackITS->GetNumClust();
950             if(trackITS->GetNumClust()) 
951                 chi2Now /= (Double_t)trackITS->GetNumClust();
952             chi2Ref = fresult->GetChi2();
953             if(fresult->GetNumClust()) 
954                 chi2Ref /= (Double_t)fresult->GetNumClust();
955             //cout<<" chi2Now and chi2Ref = "<<chi2Now<<" "<<chi2Ref<<"\n";
956             if( numClustNow > numClustRef ) {*fresult = *trackITS;} 
957             if((numClustNow == numClustRef )&& 
958                (chi2Now < chi2Ref))  {
959                 *fresult = *trackITS;
960             } // end if
961             continue;   
962         } // end if
963
964         if(trackITS->Getfnoclust()>=2)  continue;      
965         Float_t numClustNow = trackITS->GetNumClust();
966         if(numClustNow) { 
967             chi2Now = trackITS->GetChi2();
968
969             if(numClustNow<numClustRef && chi2Now>fresult->GetChi2()) continue;
970             //cout<<" chi2Now =  "<<chi2Now<<"\n";   
971               
972             chi2Now/=numClustNow;
973             if(fPtref > 1.0 && chi2Now > 30.) continue; 
974             if((fPtref >= 0.6 && fPtref<=1.0) && chi2Now > 40.) continue;
975             // if((fPtref <= 0.6 && fPtref>0.2)&& chi2Now > 40.) continue;
976             // if(fPtref <= 0.2 && chi2Now > 8.) continue;
977             if((fPtref <= 0.6 && fPtref>0.2)&& chi2Now > 30.) continue;
978             if(fPtref <= 0.2 && chi2Now > 7.) continue;     
979             /////////////////////////////
980         } // end if
981
982         Int_t layerInit = (*trackITS).GetLayer();
983         Int_t layernew = layerInit - 2;// -1 for new layer, -1 for matrix index
984         TList listoftrack;
985         Int_t ladp, ladm, detp,detm,ladinters,detinters;        
986         Int_t layerfin=layerInit-1;
987         // cout<<"Prima di intersection \n";
988         Int_t  outinters=Intersection(*trackITS,layerfin,ladinters,detinters);
989         // cout<<" outinters = "<<outinters<<"\n";
990         //  cout<<" Layer ladder detector intersection ="
991         //      <<layerfin<<" "<<ladinters<<" "<<detinters<<"\n";
992         //  cout << " phiinters zinters = "<<(*trackITS)(0) 
993         //       << " "<<(*trackITS)(1)<<"\n"; getchar();
994         if(outinters==-1) continue;
995         Int_t flaghit=0;
996         (*trackITS).SetLayer(layerfin);  // oggi
997         (*trackITS).Setfcor();           // oggi
998         if(outinters==0){
999             TVector toucLad(9), toucDet(9);      
1000             Int_t lycur=layerfin;
1001             ladp=ladinters+1;
1002             ladm=ladinters-1;
1003             if(ladm <= 0) ladm=fNlad[layerfin-1];    
1004             if(ladp > fNlad[layerfin-1]) ladp=1;  
1005             detp=detinters+1;
1006             detm=detinters-1;
1007             Int_t idetot=1;
1008             /*
1009               toucLad(0)=ladinters; toucLad(1)=ladm; toucLad(2)=ladp;
1010               toucLad(3)=ladinters; toucLad(4)=ladm; toucLad(5)=ladp;
1011               toucLad(6)=ladinters; toucLad(7)=ladm; toucLad(8)=ladp;
1012               toucDet(0)=detinters; toucDet(1)=detinters; toucDet(2)=detinters;
1013               if(detm > 0 && detp <= fNdet[layerfin-1]) {     
1014                    idetot=9;
1015                    toucDet(3)=detm; toucDet(4)=detm; toucDet(5)=detm;      
1016                    toucDet(6)=detp; toucDet(7)=detp; toucDet(8)=detp;
1017               } // end if
1018               if(detm > 0 && detp > fNdet[layerfin-1]) {   
1019                    idetot=6;
1020                    toucDet(3)=detm; toucDet(4)=detm; toucDet(5)=detm;
1021              } // end if
1022              if(detm <= 0 && detp <= fNdet[layerfin-1]) {   
1023                   idetot=6;
1024                   toucDet(3)=detp; toucDet(4)=detp; toucDet(5)=detp;
1025              } // end if
1026             */
1027             Float_t epsphi=5.0, epsz=5.0;                  
1028             if(fPtref<0.2) {epsphi=3.; epsz=3.;}     
1029             // new definition of idetot e toucLad e toucDet to be
1030             // transformed in a method
1031             // these values could be modified
1032             Float_t pigre=TMath::Pi();
1033             Float_t rangephi=5., rangez=5.;
1034             if(layerfin==1 || layerfin ==2){
1035                 rangephi=40.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1036                                                 (*trackITS).GetSigmaphi());
1037                 rangez = 40.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1038                                               (*trackITS).GetSigmaZ());
1039             } // end if
1040             if(layerfin==3 || layerfin ==4){
1041                 //rangephi=30.*fepsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1042                 //                                 (*trackITS).GetSigmaphi());
1043                 //rangez = 40.*fepsz*TMath::Sqrt(sigmazl[layerfin-1]+
1044                 //                               (*trackITS).GetSigmaZ());
1045                 rangephi=40.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1046                                                 (*trackITS).GetSigmaphi());
1047                 rangez = 50.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1048                                               (*trackITS).GetSigmaZ());
1049             } // end if
1050             if(layerfin==5 || layerfin ==6){
1051                 rangephi=20.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1052                                                 (*trackITS).GetSigmaphi());
1053                 rangez =5.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1054                                             (*trackITS).GetSigmaZ());
1055             } // end if
1056             Float_t phinters, zinters;
1057             phinters=(*trackITS).Getphi();
1058             zinters=(*trackITS).GetZ();
1059             Float_t distz = 0.0;
1060             Float_t phicm, phicp, distphim, distphip;
1061             phicm=phinters;
1062                  if(phinters>fphimax[layerfin-1][ladm-1]) phicm=phinters-2*pigre;  //corretto il 20-11-2001
1063                  distphim=TMath::Abs(phicm-fphimax[layerfin-1][ladm-1]);  //corretto il 20-11-2001
1064             phicp=phinters;
1065                  //cout<<" fNlad[layerfin-1] e ladp = "<<fNlad[layerfin-1]<<" "<<ladp<<endl;
1066                  if(phinters>fphimin[layerfin-1][ladp-1]) phicp=phinters-2.*pigre;   //corretto il 20-11-2001
1067                  distphip=TMath::Abs(phicp-fphimin[layerfin-1][ladp-1]);      //corretto il 20-11-2001
1068             Int_t flagzmin=0;
1069             Int_t flagzmax=0;
1070             idetot=1;
1071             toucLad(0)=ladinters; toucDet(0)=detinters;
1072             if(detm>0) distz=TMath::Abs(zinters-fzmax[layerfin-1][detm-1]);
1073             if(detm>0 && rangez>=distz){
1074                 flagzmin=1; 
1075                 idetot++; toucLad(idetot-1)=ladinters; toucDet(idetot-1)=detm;
1076                 if(rangephi>=distphim){
1077                     idetot++; 
1078                     toucLad(idetot-1)=ladm;
1079                     toucDet(idetot-1)=detinters;
1080                     idetot++;
1081                     toucLad(idetot-1)=ladm;
1082                     toucDet(idetot-1)=detm;
1083                 } // end if
1084                 if(rangephi>=distphip){
1085                     idetot++;
1086                     toucLad(idetot-1)=ladp;
1087                     toucDet(idetot-1)=detinters;
1088                     idetot++;
1089                     toucLad(idetot-1)=ladp;
1090                     toucDet(idetot-1)=detm;
1091                 } // end if
1092             }  //end detm>0....          
1093             if(detp<=fNdet[layerfin-1]) 
1094                 distz=TMath::Abs(zinters-fzmin[layerfin-1][detp-1]);
1095             if(detp<=fNdet[layerfin-1] && rangez>=distz){
1096                 flagzmax=1;
1097                 idetot++; toucLad(idetot-1)=ladinters; toucDet(idetot-1)=detp;
1098                 if(rangephi>=distphim){
1099                     idetot++; toucLad(idetot-1)=ladm; toucDet(idetot-1)=detp;
1100                     if(flagzmin == 0) {
1101                         idetot++;
1102                         toucLad(idetot-1)=ladm;
1103                         toucDet(idetot-1)=detinters;
1104                     } // end if
1105                 } // end if
1106                 if(rangephi>=distphip){
1107                     idetot++;
1108                     toucLad(idetot-1)=ladp;
1109                     toucDet(idetot-1)=detp;
1110                     if(flagzmin == 0) {
1111                         idetot++;
1112                         toucLad(idetot-1)=ladp; 
1113                         toucDet(idetot-1)=detinters;
1114                     } // end if
1115                 } // end if
1116             }  //end detm<fNdet[.......
1117
1118             if(flagzmin == 0 && flagzmax==0){
1119                 if(rangephi>=distphim){
1120                     idetot++; 
1121                     toucLad(idetot-1)=ladm;
1122                     toucDet(idetot-1)=detinters;
1123                 } // end if     
1124                 if(rangephi>=distphip){
1125                     idetot++; 
1126                     toucLad(idetot-1)=ladp;
1127                     toucDet(idetot-1)=detinters;
1128                 } // end if       
1129             } // end if
1130             ////////////////////////////////////////////////////////////
1131             Int_t iriv;
1132             for (iriv=0; iriv<idetot; iriv++) {  //for on detectors
1133                 ///////////////////////////////////////////////////////
1134                 /*** Rec points sorted by module *****/
1135                 /**************************************/
1136                 Int_t indexmod;       
1137                 indexmod = g1->GetModuleIndex(lycur,(Int_t)toucLad(iriv),
1138                                               (Int_t)toucDet(iriv)); 
1139                 fITS->ResetRecPoints();   
1140                 gAlice->TreeR()->GetEvent(indexmod); 
1141                 Int_t npoints=frecPoints->GetEntries();
1142        
1143                 Int_t indnew;
1144                 for(indnew=0; indnew<npoints; indnew++){
1145                     if (*(fvettid[indexmod]+indnew)==0)
1146                         recp =(AliITSRecPoint*)frecPoints->UncheckedAt(indnew);
1147                     else
1148                         continue;
1149                     TVector cluster(3),vecclust(9);
1150                     //vecclust(6)=vecclust(7)=vecclust(8)=-1.;
1151                     Double_t sigma[2];
1152   // now vecclust is with cylindrical cohordinates
1153                vecclust(0)=(Float_t)fRecCylR[indexmod][indnew];     
1154                vecclust(1)=(Float_t)fRecCylPhi[indexmod][indnew];
1155                vecclust(2)=(Float_t)fRecCylZ[indexmod][indnew];                          
1156                     vecclust(3) = (Float_t)recp->fTracks[0]; 
1157                     vecclust(4) = (Float_t)indnew;                       
1158                     vecclust(5) = (Float_t)indexmod;    
1159                     vecclust(6) = (Float_t)recp->fTracks[0];
1160                     vecclust(7) = (Float_t)recp->fTracks[1];
1161                     vecclust(8) = (Float_t)recp->fTracks[2];
1162                     sigma[0] = (Double_t)  recp->GetSigmaX2();     
1163                     sigma[1] = (Double_t) recp->GetSigmaZ2();
1164                          
1165                          cluster(0)=fRecCylR[indexmod][indnew];
1166           cluster(1)=fRecCylPhi[indexmod][indnew];
1167                          cluster(2)=fRecCylZ[indexmod][indnew];
1168                          
1169                     // cout<<" layer = "<<play<<"\n";
1170                     // cout<<" cluster prima = "<<vecclust(0)<<" "
1171                     //     <<vecclust(1)<<" "
1172                     //     <<vecclust(2)<<"\n"; getchar();    
1173
1174                     Float_t sigmatotphi, sigmatotz;  
1175                     // Float_t epsphi=5.0, epsz=5.0;                 
1176                     //if(fPtref<0.2) {epsphi=3.; epsz=3.;}
1177                     Double_t rTrack=(*trackITS).Getrtrack();
1178                     Double_t sigmaphi=sigma[0]/(rTrack*rTrack);
1179                     sigmatotphi=epsphi*TMath::Sqrt(sigmaphi + 
1180                                                    (*trackITS).GetSigmaphi());
1181                     sigmatotz=epsz*TMath::Sqrt(sigma[1] + 
1182                                                (*trackITS).GetSigmaZ());
1183                     //cout<<"cluster e sigmatotphi e track = "<<cluster(0)
1184                     //    <<" "<<cluster(1)<<" "<<sigmatotphi<<" "
1185                     //    <<vecclust(3)<<"\n";
1186                     //if(vecclust(3)==481) getchar();
1187                     if(cluster(1)<6. && (*trackITS).Getphi()>6.) 
1188                         cluster(1)=cluster(1)+(2.*TMath::Pi());
1189                     if(cluster(1)>6. && (*trackITS).Getphi()<6.) 
1190                         cluster(1)=cluster(1)-(2.*TMath::Pi());
1191                     if(TMath::Abs(cluster(1)-(*trackITS).Getphi())>sigmatotphi)
1192                         continue;
1193                     // cout<<" supero sigmaphi \n";      
1194                     AliITSTrackV1 *newTrack = new AliITSTrackV1((*trackITS));
1195                     //(*newTrack).SetLayer((*trackITS).GetLayer()-1);
1196                     if (TMath::Abs(rTrack-cluster(0))/rTrack>1e-6) 
1197                         (*newTrack).Correct(Double_t(cluster(0)));      
1198                     //cout<<" cluster(2) e(*newTrack).GetZ()="<<cluster(2)<<" "
1199                     //    << (*newTrack).GetZ()<<"\n";
1200                     if(TMath::Abs(cluster(2)-(*newTrack).GetZ()) > sigmatotz){ 
1201                         delete newTrack;
1202                         continue;
1203                     } // end if
1204                     Double_t sigmanew[2];
1205                     sigmanew[0]= sigmaphi;
1206                     sigmanew[1]=sigma[1];
1207                     Double_t m[2];
1208                     m[0]=cluster(1);
1209                     m[1]=cluster(2);
1210                     //  Double_t chi2pred=newTrack->GetPredChi2(m,sigmanew);   
1211                     // cout<<" chi2pred = "<<chi2pred<<"\n";
1212                     // if(chi2pred>fChi2max) continue; //aggiunto il 30-7-2001
1213                     if(iriv == 0) flaghit=1;
1214                     (*newTrack).AddMS();  // add the multiple scattering 
1215                                              //matrix to the covariance matrix 
1216                     (*newTrack).AddEL(1.,0);
1217
1218                     if(fflagvert){
1219                         KalmanFilterVert(newTrack,cluster,sigmanew);
1220                         //KalmanFilterVert(newTrack,cluster,sigmanew,chi2pred);
1221                     }else{
1222                         KalmanFilter(newTrack,cluster,sigmanew);
1223                     } // end if
1224                     (*newTrack).PutCluster(layernew, vecclust);
1225                     newTrack->AddClustInTrack();
1226                     listoftrack.AddLast(newTrack);
1227                 }   // end for indnew
1228             }  // end of for on detectors (iriv)
1229         }//end if(outinters==0)
1230
1231         if(flaghit==0 || outinters==-2) {
1232             AliITSTrackV1 *newTrack = new AliITSTrackV1(*trackITS);
1233             (*newTrack).Setfnoclust();           
1234             //(*newTrack).SetLayer((*trackITS).GetLayer()-1); 
1235             (*newTrack).AddMS();  // add the multiple scattering matrix
1236                                      // to the covariance matrix  
1237             (*newTrack).AddEL(1.,0);
1238             listoftrack.AddLast(newTrack);
1239         } // end if
1240
1241         //gObjectTable->Print();   // stampa memoria
1242          
1243         RecursiveTracking(&listoftrack);          
1244         listoftrack.Delete();
1245     } // end of for on tracks (index)
1246
1247     //gObjectTable->Print();   // stampa memoria
1248 }
1249 //______________________________________________________________________
1250 Int_t AliITSTrackerV1::Intersection(AliITSTrackV1 &track,Int_t layer,
1251                                     Int_t &ladder,Int_t &detector) { 
1252     // Origin  A. Badala' and G.S. Pappalardo
1253     // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1254     // Found the intersection and the detector 
1255
1256     Double_t rk=fAvrad[layer-1];
1257     if(track.DoNotCross(rk)){ /*cout<< " Do not cross \n";*/ return -1;} 
1258     track.Propagation(rk);
1259     Double_t zinters=track.GetZ();
1260     Double_t phinters=track.Getphi();
1261     //cout<<"zinters = "<<zinters<<"  phinters = "<<phinters<<"\n";
1262
1263     TVector det(9);
1264     TVector listDet(2);
1265     TVector distZCenter(2);
1266
1267     Int_t iz=0; 
1268     Int_t iD;
1269     for(iD = 1; iD<= fNdet[layer-1]; iD++) {
1270         if(zinters > fzmin[layer-1][iD-1] && zinters <= fzmax[layer-1][iD-1]) {
1271             if(iz>1) {
1272                 cout<< " Errore su iz in Intersection \n";
1273                 getchar();
1274             }else {
1275                 listDet(iz)= iD; distZCenter(iz)=TMath::Abs(zinters-det(2));
1276                 iz++;
1277             } // end if
1278         } // end if
1279     } // end for iD
1280
1281     if(iz==0) {/* cout<< " No detector along Z \n";*/ return -2;}
1282     detector=Int_t (listDet(0));
1283     if(iz>1 && (distZCenter(0)>distZCenter(1)))   detector=Int_t (listDet(1));
1284
1285     TVector listLad(2);
1286     TVector distPhiCenter(2);
1287     Int_t ip=0;
1288     Double_t pigre=TMath::Pi();
1289     Int_t iLd;   
1290     for(iLd = 1; iLd<= fNlad[layer-1]; iLd++) {
1291         Double_t phimin=fphimin[layer-1][iLd-1];
1292         Double_t phimax=fphimax[layer-1][iLd-1];
1293         Double_t phidet=fphidet[layer-1][iLd-1];
1294         Double_t phiconfr=phinters;
1295         if(phimin>phimax) {  
1296             //if(phimin <5.5) {cout<<" Error in Intersection for phi \n";
1297             // getchar();}
1298             phimin-=(2.*pigre);
1299             if(phinters>(1.5*pigre)) phiconfr=phinters-(2.*pigre); 
1300             if(phidet>(1.5*pigre)) phidet-=(2.*pigre);
1301         } // end if
1302         if(phiconfr>phimin && phiconfr<= phimax) {
1303             if(ip>1) {
1304                 cout<< " Errore su ip in Intersection \n"; getchar();
1305             }else  {
1306                 listLad(ip)= iLd;
1307                 distPhiCenter(ip)=TMath::Abs(phiconfr-phidet); ip++;
1308             } // end if
1309         } // end if
1310     } // end for iLd
1311     if(ip==0) { cout<< " No detector along phi \n"; getchar();}
1312     ladder=Int_t (listLad(0));
1313     if(ip>1 && (distPhiCenter(0)>distPhiCenter(1)))  ladder=Int_t (listLad(1));
1314     return 0;
1315 }
1316 //______________________________________________________________________
1317 void AliITSTrackerV1::KalmanFilter(AliITSTrackV1 *newTrack,TVector &cluster,
1318                                    Double_t sigma[2]){ 
1319     //Origin  A. Badala' and G.S. Pappalardo:
1320     // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1321     // Kalman filter without vertex constraint
1322     ////// Evaluation of the measurement vector ////////////////////////
1323     Double_t m[2];
1324     Double_t rk,phik,zk;
1325     rk=cluster(0);   phik=cluster(1);  zk=cluster(2);
1326     m[0]=phik;    m[1]=zk;
1327     //////////////////////// Evaluation of the error matrix V  /////////
1328     Double_t v00=sigma[0];
1329     Double_t v11=sigma[1];
1330     ////////////////////////////////////////////////////////////////////  
1331     Double_t cin00,cin10,cin20,cin30,cin40,cin11,cin21,cin31,cin41,cin22,
1332              cin32,cin42,cin33,cin43,cin44;
1333
1334     newTrack->GetCElements(cin00,
1335                            cin10,cin11,
1336                            cin20,cin21,cin22,
1337                            cin30,cin31,cin32,cin33,
1338                            cin40,cin41,cin42,cin43,cin44); //get C matrix
1339     Double_t rold00=cin00+v00;
1340     Double_t rold10=cin10;
1341     Double_t rold11=cin11+v11;
1342     ////////////////////// R matrix inversion  /////////////////////////
1343     Double_t det=rold00*rold11-rold10*rold10;
1344     Double_t r00=rold11/det;
1345     Double_t r10=-rold10/det;
1346     Double_t r11=rold00/det;
1347     ////////////////////////////////////////////////////////////////////
1348     Double_t k00=cin00*r00+cin10*r10;
1349     Double_t k01=cin00*r10+cin10*r11;
1350     Double_t k10=cin10*r00+cin11*r10;  
1351     Double_t k11=cin10*r10+cin11*r11;
1352     Double_t k20=cin20*r00+cin21*r10;  
1353     Double_t k21=cin20*r10+cin21*r11;  
1354     Double_t k30=cin30*r00+cin31*r10;  
1355     Double_t k31=cin30*r10+cin31*r11;  
1356     Double_t k40=cin40*r00+cin41*r10;
1357     Double_t k41=cin40*r10+cin41*r11;
1358     Double_t x0,x1,x2,x3,x4;
1359     newTrack->GetXElements(x0,x1,x2,x3,x4);     // get the state vector
1360     Double_t savex0=x0, savex1=x1;
1361     x0+=k00*(m[0]-savex0)+k01*(m[1]-savex1);
1362     x1+=k10*(m[0]-savex0)+k11*(m[1]-savex1);
1363     x2+=k20*(m[0]-savex0)+k21*(m[1]-savex1);
1364     x3+=k30*(m[0]-savex0)+k31*(m[1]-savex1);
1365     x4+=k40*(m[0]-savex0)+k41*(m[1]-savex1);
1366     Double_t c00,c10,c20,c30,c40,c11,c21,c31,c41,c22,c32,c42,c33,c43,c44;
1367     c00=cin00-k00*cin00-k01*cin10;
1368     c10=cin10-k00*cin10-k01*cin11;
1369     c20=cin20-k00*cin20-k01*cin21;
1370     c30=cin30-k00*cin30-k01*cin31;
1371     c40=cin40-k00*cin40-k01*cin41;
1372     c11=cin11-k10*cin10-k11*cin11;
1373     c21=cin21-k10*cin20-k11*cin21;
1374     c31=cin31-k10*cin30-k11*cin31;
1375     c41=cin41-k10*cin40-k11*cin41;
1376     c22=cin22-k20*cin20-k21*cin21;
1377     c32=cin32-k20*cin30-k21*cin31;
1378     c42=cin42-k20*cin40-k21*cin41;
1379     c33=cin33-k30*cin30-k31*cin31;
1380     c43=cin43-k30*cin40-k31*cin41;
1381     c44=cin44-k40*cin40-k41*cin41;
1382     newTrack->PutXElements(x0,x1,x2,x3,x4);  // put the new state vector
1383     newTrack->PutCElements(c00,
1384                            c10,c11,
1385                            c20,c21,c22,
1386                            c30,c31,c32,c33,
1387                            c40,c41,c42,c43,c44); // put in track the
1388                                                  // new cov matrix 
1389     Double_t vmcold00=v00-c00;
1390     Double_t vmcold10=-c10;
1391     Double_t vmcold11=v11-c11;
1392     ////////////////////// Matrix vmc inversion  ///////////////////////
1393     det=vmcold00*vmcold11-vmcold10*vmcold10;
1394     Double_t vmc00=vmcold11/det;
1395     Double_t vmc10=-vmcold10/det;
1396     Double_t vmc11=vmcold00/det;
1397     ////////////////////////////////////////////////////////////////////
1398   Double_t chi2=(m[0]-x0)*( vmc00*(m[0]-x0) + 2.*vmc10*(m[1]-x1) ) +
1399                 (m[1]-x1)*vmc11*(m[1]-x1);
1400   newTrack->SetChi2(newTrack->GetChi2()+chi2);
1401 }
1402 //----------------------------------------------------------------------
1403 //void AliITSTrackerV1::KalmanFilterVert(AliITSTrackV1 *newTrack,
1404 //                                       TVector &cluster,Double_t sigma[2]){
1405 void AliITSTrackerV1::KalmanFilterVert(AliITSTrackV1 *newTrack,
1406                                        TVector &cluster,Double_t sigma[2]
1407                                        /*, Double_t chi2pred*/){
1408     //Origin  A. Badala' and G.S. Pappalardo:
1409     // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it 
1410     // Kalman filter with vertex constraint
1411     ///////////////////// Evaluation of the measurement vector m ////////
1412     Double_t m[4];
1413     Double_t rk,phik,zk;
1414     rk=cluster(0);   phik=cluster(1);  zk=cluster(2);
1415     m[0]=phik;    m[1]=zk;
1416     Double_t cc=(*newTrack).GetC();
1417     Double_t zv=(*newTrack).GetZv(); 
1418     Double_t dv=(*newTrack).GetDv();
1419     Double_t cy=cc/2.;
1420     Double_t tgl= (zk-zv)*cy/TMath::ASin(cy*rk);
1421     m[2]=dv;    m[3]=tgl;
1422     /////////////////////// Evaluation of the error matrix V  //////////
1423     Int_t layer=newTrack->GetLayer();
1424     Double_t v00=sigma[0];
1425     Double_t v11=sigma[1];
1426     Double_t v31=sigma[1]/rk;
1427     Double_t sigmaDv=newTrack->GetsigmaDv();
1428     Double_t v22=sigmaDv*sigmaDv  + newTrack->Getd2(layer-1);
1429     Double_t v32=newTrack->Getdtgl(layer-1);
1430     Double_t sigmaZv=newTrack->GetsigmaZv();  
1431     Double_t v33=(sigma[1]+sigmaZv*sigmaZv)/(rk*rk)+newTrack->Gettgl2(layer-1);
1432     //////////////////////////////////////////////////////////////////
1433     Double_t cin00,cin10,cin11,cin20,cin21,cin22,
1434              cin30,cin31,cin32,cin33,cin40,cin41,cin42,cin43,cin44;
1435     newTrack->GetCElements(cin00,
1436                            cin10,cin11,
1437                            cin20,cin21,cin22,
1438                            cin30,cin31,cin32,cin33,
1439                            cin40,cin41,cin42,cin43,cin44); //get C matrix
1440     Double_t r[4][4];
1441     r[0][0]=cin00+v00;
1442     r[1][0]=cin10;
1443     r[2][0]=cin20;
1444     r[3][0]=cin30;
1445     r[1][1]=cin11+v11;
1446     r[2][1]=cin21;
1447     r[3][1]=cin31+sigma[1]/rk;
1448     r[2][2]=cin22+sigmaDv*sigmaDv+newTrack->Getd2(layer-1);
1449     r[3][2]=cin32+newTrack->Getdtgl(layer-1);
1450     r[3][3]=cin33+(sigma[1]+sigmaZv*sigmaZv)/(rk*rk)+
1451                                     newTrack->Gettgl2(layer-1);
1452     r[0][1]=r[1][0]; r[0][2]=r[2][0]; r[0][3]=r[3][0]; 
1453     r[1][2]=r[2][1]; r[1][3]=r[3][1]; r[2][3]=r[3][2];
1454     /////////////////////  Matrix R inversion //////////////////////////
1455     const Int_t kn=4;
1456     Double_t big, hold;
1457     Double_t d=1.;
1458     Int_t ll[kn],mm[kn];
1459     Int_t i,j,k;
1460
1461     for(k=0; k<kn; k++) {
1462         ll[k]=k;
1463         mm[k]=k;
1464         big=r[k][k];
1465         for(j=k; j<kn ; j++) {
1466             for (i=j; i<kn; i++) {
1467                 if(TMath::Abs(big) < TMath::Abs(r[i][j]) ) { 
1468                     big=r[i][j]; 
1469                     ll[k]=i; 
1470                     mm[k]=j; 
1471                 } // end if
1472             } // end for i
1473         } // end for j
1474         //
1475         j= ll[k];
1476         if(j > k) {
1477             for(i=0; i<kn; i++) {
1478                 hold=-r[k][i]; 
1479                 r[k][i]=r[j][i]; 
1480                 r[j][i]=hold;
1481             } // end for i
1482         }// end if
1483         //
1484         i=mm[k];
1485         if(i > k ) { 
1486             for(j=0; j<kn; j++) {
1487                 hold=-r[j][k];
1488                 r[j][k]=r[j][i]; 
1489                 r[j][i]=hold; 
1490             } // end for j
1491         } // end if
1492         //
1493         if(!big) {
1494             d=0.;
1495             cout << "Singular matrix\n"; 
1496         } // end if
1497         for(i=0; i<kn; i++) {
1498             if(i == k) { continue; }    
1499             r[i][k]=r[i][k]/(-big);
1500         } // end for i
1501         //
1502         for(i=0; i<kn; i++) {
1503             hold=r[i][k];
1504             for(j=0; j<kn; j++) {
1505                 if(i == k || j == k) continue;
1506                 r[i][j]=hold*r[k][j]+r[i][j];
1507             } // end for j
1508         } // end of ri
1509         // 
1510         for(j=0; j<kn; j++) {
1511             if(j == k) continue;
1512             r[k][j]=r[k][j]/big;
1513         } // end for j
1514         //
1515         d=d*big;
1516         //
1517         r[k][k]=1./big;        
1518     } // end for k
1519     //
1520     for(k=kn-1; k>=0; k--) {
1521         i=ll[k];
1522         if(i > k) {
1523             for (j=0; j<kn; j++) {
1524                 hold=r[j][k];
1525                 r[j][k]=-r[j][i]; 
1526                 r[j][i]=hold;
1527             } // end for j
1528         } // end if i
1529         j=mm[k];
1530         if(j > k) {
1531             for (i=0; i<kn; i++) {
1532                 hold=r[k][i]; 
1533                 r[k][i]=-r[j][i]; 
1534                 r[j][i]=hold;
1535             } // end for i
1536         } // end if
1537     } // end for k
1538     ////////////////////////////////////////////////////////////////////
1539     Double_t k00=cin00*r[0][0]+cin10*r[1][0]+cin20*r[2][0]+cin30*r[3][0];
1540     Double_t k01=cin00*r[1][0]+cin10*r[1][1]+cin20*r[2][1]+cin30*r[3][1];
1541     Double_t k02=cin00*r[2][0]+cin10*r[2][1]+cin20*r[2][2]+cin30*r[3][2];
1542     Double_t k03=cin00*r[3][0]+cin10*r[3][1]+cin20*r[3][2]+cin30*r[3][3];
1543     Double_t k10=cin10*r[0][0]+cin11*r[1][0]+cin21*r[2][0]+cin31*r[3][0];  
1544     Double_t k11=cin10*r[1][0]+cin11*r[1][1]+cin21*r[2][1]+cin31*r[3][1];
1545     Double_t k12=cin10*r[2][0]+cin11*r[2][1]+cin21*r[2][2]+cin31*r[3][2];
1546     Double_t k13=cin10*r[3][0]+cin11*r[3][1]+cin21*r[3][2]+cin31*r[3][3];
1547     Double_t k20=cin20*r[0][0]+cin21*r[1][0]+cin22*r[2][0]+cin32*r[3][0];  
1548     Double_t k21=cin20*r[1][0]+cin21*r[1][1]+cin22*r[2][1]+cin32*r[3][1];  
1549     Double_t k22=cin20*r[2][0]+cin21*r[2][1]+cin22*r[2][2]+cin32*r[3][2];
1550     Double_t k23=cin20*r[3][0]+cin21*r[3][1]+cin22*r[3][2]+cin32*r[3][3];
1551     Double_t k30=cin30*r[0][0]+cin31*r[1][0]+cin32*r[2][0]+cin33*r[3][0];  
1552     Double_t k31=cin30*r[1][0]+cin31*r[1][1]+cin32*r[2][1]+cin33*r[3][1];  
1553     Double_t k32=cin30*r[2][0]+cin31*r[2][1]+cin32*r[2][2]+cin33*r[3][2];  
1554     Double_t k33=cin30*r[3][0]+cin31*r[3][1]+cin32*r[3][2]+cin33*r[3][3];
1555     Double_t k40=cin40*r[0][0]+cin41*r[1][0]+cin42*r[2][0]+cin43*r[3][0];
1556     Double_t k41=cin40*r[1][0]+cin41*r[1][1]+cin42*r[2][1]+cin43*r[3][1];
1557     Double_t k42=cin40*r[2][0]+cin41*r[2][1]+cin42*r[2][2]+cin43*r[3][2];  
1558     Double_t k43=cin40*r[3][0]+cin41*r[3][1]+cin42*r[3][2]+cin43*r[3][3];
1559
1560     Double_t x0,x1,x2,x3,x4;
1561     newTrack->GetXElements(x0,x1,x2,x3,x4);     // get the state vector
1562     Double_t savex0=x0, savex1=x1, savex2=x2, savex3=x3;
1563     x0+=k00*(m[0]-savex0)+k01*(m[1]-savex1)+k02*(m[2]-savex2)+
1564         k03*(m[3]-savex3);
1565     x1+=k10*(m[0]-savex0)+k11*(m[1]-savex1)+k12*(m[2]-savex2)+
1566         k13*(m[3]-savex3);
1567     x2+=k20*(m[0]-savex0)+k21*(m[1]-savex1)+k22*(m[2]-savex2)+
1568         k23*(m[3]-savex3);
1569     x3+=k30*(m[0]-savex0)+k31*(m[1]-savex1)+k32*(m[2]-savex2)+
1570         k33*(m[3]-savex3);
1571     x4+=k40*(m[0]-savex0)+k41*(m[1]-savex1)+k42*(m[2]-savex2)+
1572         k43*(m[3]-savex3);
1573     Double_t c00,c10,c20,c30,c40,c11,c21,c31,c41,c22,c32,c42,c33,c43,c44;
1574     c00=cin00-k00*cin00-k01*cin10-k02*cin20-k03*cin30;
1575     c10=cin10-k00*cin10-k01*cin11-k02*cin21-k03*cin31;
1576     c20=cin20-k00*cin20-k01*cin21-k02*cin22-k03*cin32;
1577     c30=cin30-k00*cin30-k01*cin31-k02*cin32-k03*cin33;
1578     c40=cin40-k00*cin40-k01*cin41-k02*cin42-k03*cin43;
1579     c11=cin11-k10*cin10-k11*cin11-k12*cin21-k13*cin31;
1580     c21=cin21-k10*cin20-k11*cin21-k12*cin22-k13*cin32;
1581     c31=cin31-k10*cin30-k11*cin31-k12*cin32-k13*cin33;
1582     c41=cin41-k10*cin40-k11*cin41-k12*cin42-k13*cin43;
1583     c22=cin22-k20*cin20-k21*cin21-k22*cin22-k23*cin32;
1584     c32=cin32-k20*cin30-k21*cin31-k22*cin32-k23*cin33;
1585     c42=cin42-k20*cin40-k21*cin41-k22*cin42-k23*cin43;
1586     c33=cin33-k30*cin30-k31*cin31-k32*cin32-k33*cin33;
1587     c43=cin43-k30*cin40-k31*cin41-k32*cin42-k33*cin43;
1588     c44=cin44-k40*cin40-k41*cin41-k42*cin42-k43*cin43;
1589
1590     newTrack->PutXElements(x0,x1,x2,x3,x4); // put the new state vector
1591     newTrack->PutCElements(c00,
1592                            c10,c11,
1593                            c20,c21,c22,
1594                            c30,c31,c32,c33,
1595                            c40,c41,c42,c43,c44); // put in track the
1596                                                  // new cov matrix
1597     Double_t vmc[4][4];
1598     vmc[0][0]=v00-c00; vmc[1][0]=-c10; vmc[2][0]=-c20; vmc[3][0]=-c30;
1599     vmc[1][1]=v11-c11; vmc[2][1]=-c21; vmc[3][1]=v31-c31;
1600     vmc[2][2]=v22-c22; vmc[3][2]=v32-c32;
1601     vmc[3][3]=v33-c33;
1602     vmc[0][1]=vmc[1][0]; vmc[0][2]=vmc[2][0]; vmc[0][3]=vmc[3][0];
1603     vmc[1][2]=vmc[2][1]; vmc[1][3]=vmc[3][1];
1604     vmc[2][3]=vmc[3][2];
1605     /////////////////////// vmc matrix inversion /////////////////////// 
1606     d=1.;
1607     for(k=0; k<kn; k++) {
1608         ll[k]=k;
1609         mm[k]=k;
1610         big=vmc[k][k];
1611         for(j=k; j<kn ; j++) {
1612             for (i=j; i<kn; i++) {
1613                 if(TMath::Abs(big) < TMath::Abs(vmc[i][j]) ) {
1614                     big=vmc[i][j]; 
1615                     ll[k]=i;
1616                     mm[k]=j;
1617                 } // end if
1618             } // end for i
1619         } // end for j
1620         //
1621         j= ll[k];
1622         if(j > k) {
1623             for(i=0; i<kn; i++) {
1624                 hold=-vmc[k][i]; 
1625                 vmc[k][i]=vmc[j][i];
1626                 vmc[j][i]=hold;
1627             } // end for i
1628         } // end if
1629         //
1630         i=mm[k];
1631         if(i > k ) { 
1632             for(j=0; j<kn; j++) { 
1633                 hold=-vmc[j][k]; 
1634                 vmc[j][k]=vmc[j][i]; 
1635                 vmc[j][i]=hold; 
1636             } // end for j
1637         } // end if
1638         //
1639         if(!big) {
1640             d=0.;
1641             cout << "Singular matrix\n"; 
1642         } // end if
1643         for(i=0; i<kn; i++) {
1644             if(i == k) continue;    
1645             vmc[i][k]=vmc[i][k]/(-big);
1646         }    // end for i
1647         //
1648         for(i=0; i<kn; i++) {
1649             hold=vmc[i][k];
1650             for(j=0; j<kn; j++) {
1651                 if(i == k || j == k) continue;
1652                 vmc[i][j]=hold*vmc[k][j]+vmc[i][j];
1653             } // end for j
1654         } // end for i
1655         //  
1656         for(j=0; j<kn; j++) {
1657             if(j == k) continue;
1658             vmc[k][j]=vmc[k][j]/big;
1659         } // end for j
1660         //
1661         d=d*big;
1662         //
1663         vmc[k][k]=1./big;        
1664     }  // end for k
1665     //  
1666     for(k=kn-1; k>=0; k--) {
1667         i=ll[k];
1668         if(i > k) {
1669             for (j=0; j<kn; j++) {
1670                 hold=vmc[j][k]; 
1671                 vmc[j][k]=-vmc[j][i];
1672                 vmc[j][i]=hold;
1673             } // end for j
1674         } // end if i>k
1675         j=mm[k];
1676         if(j > k) {
1677             for (i=0; i<kn; i++) {
1678                 hold=vmc[k][i]; 
1679                 vmc[k][i]=-vmc[j][i]; 
1680                 vmc[j][i]=hold;
1681             } // end for i
1682         } // end if j>k
1683     } // end for k
1684     ////////////////////////////////////////////////////////////////////
1685     Double_t chi2=(m[0]-x0)*( vmc[0][0]*(m[0]-x0) + 2.*vmc[1][0]*(m[1]-x1) + 
1686                               2.*vmc[2][0]*(m[2]-x2)+ 2.*vmc[3][0]*(m[3]-x3) )+
1687         (m[1]-x1)* ( vmc[1][1]*(m[1]-x1) + 2.*vmc[2][1]*(m[2]-x2)+ 
1688                      2.*vmc[3][1]*(m[3]-x3) ) +
1689         (m[2]-x2)* ( vmc[2][2]*(m[2]-x2)+ 2.*vmc[3][2]*(m[3]-x3) ) +
1690         (m[3]-x3)*vmc[3][3]*(m[3]-x3);
1691     //cout<<" chi2 kalman = "<<chi2<<"\n";  getchar(); 
1692     newTrack->SetChi2(newTrack->GetChi2()+chi2);   
1693     //   newTrack->SetChi2(newTrack->GetChi2()+chi2pred);
1694 }
1695