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