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