]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPC.cxx
Coding convention rules obeyed
[u/mrichter/AliRoot.git] / TPC / AliTPC.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 /*
17 $Log$
18 Revision 1.24  2000/10/05 16:06:09  kowal2
19 Forward declarations. Changes due to a new class AliComplexCluster.
20
21 Revision 1.23  2000/10/02 21:28:18  fca
22 Removal of useless dependecies via forward declarations
23
24 Revision 1.22  2000/07/10 20:57:39  hristov
25 Update of TPC code and macros by M.Kowalski
26
27 Revision 1.19.2.4  2000/06/26 07:39:42  kowal2
28 Changes to obey the coding rules
29
30 Revision 1.19.2.3  2000/06/25 08:38:41  kowal2
31 Splitted from AliTPCtracking
32
33 Revision 1.19.2.2  2000/06/16 12:59:28  kowal2
34 Changed parameter settings
35
36 Revision 1.19.2.1  2000/06/09 07:15:07  kowal2
37
38 Defaults loaded automatically (hard-wired)
39 Optional parameters can be set via macro called in the constructor
40
41 Revision 1.19  2000/04/18 19:00:59  fca
42 Small bug fixes to TPC files
43
44 Revision 1.18  2000/04/17 09:37:33  kowal2
45 removed obsolete AliTPCDigitsDisplay.C
46
47 Revision 1.17.2.2  2000/04/10 08:15:12  kowal2
48
49 New, experimental data structure from M. Ivanov
50 New tracking algorithm
51 Different pad geometry for different sectors
52 Digitization rewritten
53
54 Revision 1.17.2.1  2000/04/10 07:56:53  kowal2
55 Not used anymore - removed
56
57 Revision 1.17  2000/01/19 17:17:30  fca
58 Introducing a list of lists of hits -- more hits allowed for detector now
59
60 Revision 1.16  1999/11/05 09:29:23  fca
61 Accept only signals > 0
62
63 Revision 1.15  1999/10/08 06:26:53  fca
64 Removed ClustersIndex - not used anymore
65
66 Revision 1.14  1999/09/29 09:24:33  fca
67 Introduction of the Copyright and cvs Log
68
69 */
70
71 ///////////////////////////////////////////////////////////////////////////////
72 //                                                                           //
73 //  Time Projection Chamber                                                  //
74 //  This class contains the basic functions for the Time Projection Chamber  //
75 //  detector. Functions specific to one particular geometry are              //
76 //  contained in the derived classes                                         //
77 //                                                                           //
78 //Begin_Html
79 /*
80 <img src="picts/AliTPCClass.gif">
81 */
82 //End_Html
83 //                                                                           //
84 //                                                                          //
85 ///////////////////////////////////////////////////////////////////////////////
86
87 //
88
89 #include <TMath.h>
90 #include <TRandom.h>
91 #include <TVector.h>
92 #include <TMatrix.h>
93 #include <TGeometry.h>
94 #include <TNode.h>
95 #include <TTUBS.h>
96 #include <TObjectTable.h>
97 #include "TParticle.h"
98 #include "AliTPC.h"
99 #include <TFile.h>       
100 #include "AliRun.h"
101 #include <iostream.h>
102 #include <fstream.h>
103 #include "AliMC.h"
104 #include "AliMagF.h"
105
106
107 #include "AliTPCParamSR.h"
108 #include "AliTPCPRF2D.h"
109 #include "AliTPCRF1D.h"
110 #include "AliDigits.h"
111 #include "AliSimDigits.h"
112 #include "AliTPCTrackHits.h"
113 #include "AliPoints.h"
114 #include "AliArrayBranch.h"
115
116
117 #include "AliTPCDigitsArray.h"
118 #include "AliComplexCluster.h"
119 #include "AliClusters.h"
120 #include "AliTPCClustersRow.h"
121 #include "AliTPCClustersArray.h"
122
123 #include "AliTPCcluster.h"
124 #include "AliTPCclusterer.h"
125 #include "AliTPCtracker.h"
126
127 #include <TInterpreter.h>
128 #include <TTree.h>
129
130
131
132 ClassImp(AliTPC) 
133
134 //_____________________________________________________________________________
135 AliTPC::AliTPC()
136 {
137   //
138   // Default constructor
139   //
140   fIshunt   = 0;
141   fHits     = 0;
142   fDigits   = 0;
143   fNsectors = 0;
144   //MI changes
145   fDigitsArray = 0;
146   fClustersArray = 0;
147   fTPCParam=0;
148   fTrackHits = 0;  
149   fHitType = 2;  
150   fTPCParam = 0; 
151 }
152  
153 //_____________________________________________________________________________
154 AliTPC::AliTPC(const char *name, const char *title)
155       : AliDetector(name,title)
156 {
157   //
158   // Standard constructor
159   //
160
161   //
162   // Initialise arrays of hits and digits 
163   fHits     = new TClonesArray("AliTPChit",  176);
164   gAlice->AddHitList(fHits);
165   //MI change  
166   fDigitsArray = 0;
167   fClustersArray= 0;
168   //
169   fTrackHits = new AliTPCTrackHits;  //MI - 13.09.2000
170   fTrackHits->SetHitPrecision(0.002);
171   fTrackHits->SetStepPrecision(0.003);  
172   fTrackHits->SetMaxDistance(100); 
173
174   fHitType = 2;
175   //
176   // Initialise counters
177   fNsectors = 0;
178
179   //
180   fIshunt     =  0;
181   //
182   // Initialise color attributes
183   SetMarkerColor(kYellow);
184
185   //
186   //  Set TPC parameters
187   //
188
189   if (!strcmp(title,"Default")) {  
190      fTPCParam = new AliTPCParamSR;
191   } else {
192     cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
193     fTPCParam=0;
194   }
195
196 }
197
198 //_____________________________________________________________________________
199 AliTPC::~AliTPC()
200 {
201   //
202   // TPC destructor
203   //
204
205   fIshunt   = 0;
206   delete fHits;
207   delete fDigits;
208   delete fTPCParam;
209   delete fTrackHits; //MI 15.09.2000
210 }
211
212 //_____________________________________________________________________________
213 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
214 {
215   //
216   // Add a hit to the list
217   //
218   //  TClonesArray &lhits = *fHits;
219   //  new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
220   if (fHitType&1){
221     TClonesArray &lhits = *fHits;
222     new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
223   }
224   if (fHitType&2)
225    AddHit2(track,vol,hits);
226 }
227  
228 //_____________________________________________________________________________
229 void AliTPC::BuildGeometry()
230 {
231
232   //
233   // Build TPC ROOT TNode geometry for the event display
234   //
235   TNode *nNode, *nTop;
236   TTUBS *tubs;
237   Int_t i;
238   const int kColorTPC=19;
239   char name[5], title[25];
240   const Double_t kDegrad=TMath::Pi()/180;
241   const Double_t kRaddeg=180./TMath::Pi();
242
243
244   Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
245   Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
246
247   Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
248   Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
249
250   Int_t nLo = fTPCParam->GetNInnerSector()/2;
251   Int_t nHi = fTPCParam->GetNOuterSector()/2;  
252
253   const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
254   const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
255   const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
256   const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);  
257
258
259   const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
260   const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
261
262   Double_t rl,ru;
263   
264
265   //
266   // Get ALICE top node
267   //
268
269   nTop=gAlice->GetGeometry()->GetNode("alice");
270
271   //  inner sectors
272
273   rl = fTPCParam->GetInnerRadiusLow();
274   ru = fTPCParam->GetInnerRadiusUp();
275  
276
277   for(i=0;i<nLo;i++) {
278     sprintf(name,"LS%2.2d",i);
279     name[4]='\0';
280     sprintf(title,"TPC low sector %3d",i);
281     title[24]='\0';
282     
283     tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
284                      kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
285     tubs->SetNumberOfDivisions(1);
286     nTop->cd();
287     nNode = new TNode(name,title,name,0,0,0,"");
288     nNode->SetLineColor(kColorTPC);
289     fNodes->Add(nNode);
290   }
291
292   // Outer sectors
293
294   rl = fTPCParam->GetOuterRadiusLow();
295   ru = fTPCParam->GetOuterRadiusUp();
296
297   for(i=0;i<nHi;i++) {
298     sprintf(name,"US%2.2d",i);
299     name[4]='\0';
300     sprintf(title,"TPC upper sector %d",i);
301     title[24]='\0';
302     tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
303                      khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
304     tubs->SetNumberOfDivisions(1);
305     nTop->cd();
306     nNode = new TNode(name,title,name,0,0,0,"");
307     nNode->SetLineColor(kColorTPC);
308     fNodes->Add(nNode);
309   }
310
311 }    
312
313 //_____________________________________________________________________________
314 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
315 {
316   //
317   // Calculate distance from TPC to mouse on the display
318   // Dummy procedure
319   //
320   return 9999;
321 }
322
323 void AliTPC::Clusters2Tracks(TFile *of) {
324   //-----------------------------------------------------------------
325   // This is a track finder.
326   //-----------------------------------------------------------------
327   AliTPCtracker::Clusters2Tracks(fTPCParam,of);
328 }
329
330 //_____________________________________________________________________________
331 void AliTPC::CreateMaterials()
332 {
333   //-----------------------------------------------
334   // Create Materials for for TPC simulations
335   //-----------------------------------------------
336
337   //-----------------------------------------------------------------
338   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
339   //-----------------------------------------------------------------
340
341   Int_t iSXFLD=gAlice->Field()->Integ();
342   Float_t sXMGMX=gAlice->Field()->Max();
343
344   Float_t amat[5]; // atomic numbers
345   Float_t zmat[5]; // z
346   Float_t wmat[5]; // proportions
347
348   Float_t density;
349   Float_t apure[2];
350
351
352   //***************** Gases *************************
353   
354   //-------------------------------------------------
355   // pure gases
356   //-------------------------------------------------
357
358   // Neon
359
360
361   amat[0]= 20.18;
362   zmat[0]= 10.;  
363   density = 0.0009;
364  
365   apure[0]=amat[0];
366
367   AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
368
369   // Argon
370
371   amat[0]= 39.948;
372   zmat[0]= 18.;  
373   density = 0.001782;  
374
375   apure[1]=amat[0];
376
377   AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
378  
379
380   //--------------------------------------------------------------
381   // gases - compounds
382   //--------------------------------------------------------------
383
384   Float_t amol[3];
385
386   // CO2
387
388   amat[0]=12.011;
389   amat[1]=15.9994;
390
391   zmat[0]=6.;
392   zmat[1]=8.;
393
394   wmat[0]=1.;
395   wmat[1]=2.;
396
397   density=0.001977;
398
399   amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
400
401   AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
402   
403   // CF4
404
405   amat[0]=12.011;
406   amat[1]=18.998;
407
408   zmat[0]=6.;
409   zmat[1]=9.;
410  
411   wmat[0]=1.;
412   wmat[1]=4.;
413  
414   density=0.003034;
415
416   amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
417
418   AliMixture(11,"CF4",amat,zmat,density,-2,wmat); 
419
420
421   // CH4
422
423   amat[0]=12.011;
424   amat[1]=1.;
425
426   zmat[0]=6.;
427   zmat[1]=1.;
428
429   wmat[0]=1.;
430   wmat[1]=4.;
431
432   density=0.000717;
433
434   amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
435
436   AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
437
438   //----------------------------------------------------------------
439   // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
440   //----------------------------------------------------------------
441
442   char namate[21]; 
443   density = 0.;
444   Float_t am=0;
445   Int_t nc;
446   Float_t rho,absl,X0,buf[1];
447   Int_t nbuf;
448   Float_t a,z;
449
450   for(nc = 0;nc<fNoComp;nc++)
451     {
452     
453       // retrive material constants
454       
455       gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
456
457       amat[nc] = a;
458       zmat[nc] = z;
459
460       Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
461  
462       am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]); 
463       density += fMixtProp[nc]*rho;  // density of the mixture
464       
465     }
466
467   // mixture proportions by weight!
468
469   for(nc = 0;nc<fNoComp;nc++)
470     {
471
472       Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
473
474       wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ? 
475                  apure[nnc] : amol[nnc])/am;
476
477     } 
478
479   // Drift gases 1 - nonsensitive, 2 - sensitive
480
481   AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
482   AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
483
484
485   // Air
486
487   amat[0] = 14.61;
488   zmat[0] = 7.3;
489   density = 0.001205;
490
491   AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.); 
492
493
494   //----------------------------------------------------------------------
495   //               solid materials
496   //----------------------------------------------------------------------
497
498
499   // Kevlar C14H22O2N2
500
501   amat[0] = 12.011;
502   amat[1] = 1.;
503   amat[2] = 15.999;
504   amat[3] = 14.006;
505
506   zmat[0] = 6.;
507   zmat[1] = 1.;
508   zmat[2] = 8.;
509   zmat[3] = 7.;
510
511   wmat[0] = 14.;
512   wmat[1] = 22.;
513   wmat[2] = 2.;
514   wmat[3] = 2.;
515
516   density = 1.45;
517
518   AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);  
519
520   // NOMEX
521
522   amat[0] = 12.011;
523   amat[1] = 1.;
524   amat[2] = 15.999;
525   amat[3] = 14.006;
526
527   zmat[0] = 6.;
528   zmat[1] = 1.;
529   zmat[2] = 8.;
530   zmat[3] = 7.;
531
532   wmat[0] = 14.;
533   wmat[1] = 22.;
534   wmat[2] = 2.;
535   wmat[3] = 2.;
536
537   density = 0.03;
538
539   
540   AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
541
542   // Makrolon C16H18O3
543
544   amat[0] = 12.011;
545   amat[1] = 1.;
546   amat[2] = 15.999;
547
548   zmat[0] = 6.;
549   zmat[1] = 1.;
550   zmat[2] = 8.;
551
552   wmat[0] = 16.;
553   wmat[1] = 18.;
554   wmat[2] = 3.;
555   
556   density = 1.2;
557
558   AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
559   
560   // Mylar C5H4O2
561
562   amat[0]=12.011;
563   amat[1]=1.;
564   amat[2]=15.9994;
565
566   zmat[0]=6.;
567   zmat[1]=1.;
568   zmat[2]=8.;
569
570   wmat[0]=5.;
571   wmat[1]=4.;
572   wmat[2]=2.; 
573
574   density = 1.39;
575   
576   AliMixture(37, "Mylar",amat,zmat,density,-3,wmat); 
577
578   // G10 60% SiO2 + 40% epoxy, I use A and Z for SiO2
579
580   amat[0]=28.086;
581   amat[1]=15.9994;
582
583   zmat[0]=14.;
584   zmat[1]=8.;
585
586   wmat[0]=1.;
587   wmat[1]=2.;
588
589   density = 1.7;
590
591   AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz
592  
593   gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf); 
594
595   AliMaterial(39,"G10",amat[0],zmat[0],density,999.,999.); 
596
597   // Al
598
599   amat[0] = 26.98;
600   zmat[0] = 13.;
601
602   density = 2.7;
603
604   AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
605
606   // Si
607
608   amat[0] = 28.086;
609   zmat[0] = 14.,
610
611   density = 2.33;
612
613   AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
614
615   // Cu
616
617   amat[0] = 63.546;
618   zmat[0] = 29.;
619
620   density = 8.96;
621
622   AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
623
624   // Tedlar C2H3F
625
626   amat[0] = 12.011;
627   amat[1] = 1.;
628   amat[2] = 18.998;
629
630   zmat[0] = 6.;
631   zmat[1] = 1.;
632   zmat[2] = 9.;
633
634   wmat[0] = 2.;
635   wmat[1] = 3.; 
636   wmat[2] = 1.;
637
638   density = 1.71;
639
640   AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);  
641
642
643   // Plexiglas  C5H8O2
644
645   amat[0]=12.011;
646   amat[1]=1.;
647   amat[2]=15.9994;
648
649   zmat[0]=6.;
650   zmat[1]=1.;
651   zmat[2]=8.;
652
653   wmat[0]=5.;
654   wmat[1]=8.;
655   wmat[2]=2.;
656
657   density=1.18;
658
659   AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
660
661
662
663   //----------------------------------------------------------
664   // tracking media for gases
665   //----------------------------------------------------------
666
667   AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
668   AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
669   AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
670   AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001); 
671
672   //-----------------------------------------------------------  
673   // tracking media for solids
674   //-----------------------------------------------------------
675   
676   AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
677   AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
678   AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
679   AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
680   AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
681   AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
682   AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
683   AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
684   AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
685   AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
686     
687 }
688
689
690 void AliTPC::Digits2Clusters(TFile *of)
691 {
692   //-----------------------------------------------------------------
693   // This is a simple cluster finder.
694   //-----------------------------------------------------------------
695   AliTPCclusterer::Digits2Clusters(fTPCParam,of);
696 }
697
698 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
699 extern Double_t SigmaZ2(Double_t, Double_t);
700 //_____________________________________________________________________________
701 void AliTPC::Hits2Clusters(TFile *of)
702 {
703   //--------------------------------------------------------
704   // TPC simple cluster generator from hits
705   // obtained from the TPC Fast Simulator
706   // The point errors are taken from the parametrization
707   //--------------------------------------------------------
708
709   //-----------------------------------------------------------------
710   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
711   //-----------------------------------------------------------------
712   // Adopted to Marian's cluster data structure by I.Belikov, CERN,
713   // Jouri.Belikov@cern.ch
714   //----------------------------------------------------------------
715   
716   /////////////////////////////////////////////////////////////////////////////
717   //
718   //---------------------------------------------------------------------
719   //   ALICE TPC Cluster Parameters
720   //--------------------------------------------------------------------
721        
722   
723
724   // Cluster width in rphi
725   const Float_t kACrphi=0.18322;
726   const Float_t kBCrphi=0.59551e-3;
727   const Float_t kCCrphi=0.60952e-1;
728   // Cluster width in z
729   const Float_t kACz=0.19081;
730   const Float_t kBCz=0.55938e-3;
731   const Float_t kCCz=0.30428;
732
733   TDirectory *savedir=gDirectory; 
734
735   if (!of->IsOpen()) {
736      cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
737      return;
738   }
739
740    if(fTPCParam == 0){
741      printf("AliTPCParam MUST be created firstly\n");
742      return;
743    }
744
745   Float_t sigmaRphi,sigmaZ,clRphi,clZ;
746   //
747   TParticle *particle; // pointer to a given particle
748   AliTPChit *tpcHit; // pointer to a sigle TPC hit
749   TClonesArray *particles; //pointer to the particle list
750   Int_t sector;
751   Int_t ipart;
752   Float_t xyz[5];
753   Float_t pl,pt,tanth,rpad,ratio;
754   Float_t cph,sph;
755   
756   //---------------------------------------------------------------
757   //  Get the access to the tracks 
758   //---------------------------------------------------------------
759   
760   TTree *tH = gAlice->TreeH();
761   Stat_t ntracks = tH->GetEntries();
762   particles=gAlice->Particles();
763
764   //Switch to the output file
765   of->cd();
766
767   fTPCParam->Write(fTPCParam->GetTitle());
768   AliTPCClustersArray carray;
769   carray.Setup(fTPCParam);
770   carray.SetClusterType("AliTPCcluster");
771   carray.MakeTree();
772
773   Int_t nclusters=0; //cluster counter
774   
775   //------------------------------------------------------------
776   // Loop over all sectors (72 sectors for 20 deg
777   // segmentation for both lower and upper sectors)
778   // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
779   // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
780   //
781   // First cluster for sector 0 starts at "0"
782   //------------------------------------------------------------
783    
784   for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
785     //MI change
786     fTPCParam->AdjustCosSin(isec,cph,sph);
787     
788     //------------------------------------------------------------
789     // Loop over tracks
790     //------------------------------------------------------------
791     
792     for(Int_t track=0;track<ntracks;track++){
793       ResetHits();
794       tH->GetEvent(track);
795       //
796       //  Get number of the TPC hits
797       //
798       // nhits=fHits->GetEntriesFast();
799       //
800      
801        tpcHit = (AliTPChit*)FirstHit(-1);
802
803       // Loop over hits
804       //
805        //   for(Int_t hit=0;hit<nhits;hit++){
806        //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
807
808        while(tpcHit){
809  
810          if (tpcHit->fQ == 0.) {
811            tpcHit = (AliTPChit*) NextHit();
812            continue; //information about track (I.Belikov)
813          }
814         sector=tpcHit->fSector; // sector number
815
816
817         //      if(sector != isec) continue; //terminate iteration
818
819        if(sector != isec){
820          tpcHit = (AliTPChit*) NextHit();
821          continue; 
822        }
823         ipart=tpcHit->Track();
824         particle=(TParticle*)particles->UncheckedAt(ipart);
825         pl=particle->Pz();
826         pt=particle->Pt();
827         if(pt < 1.e-9) pt=1.e-9;
828         tanth=pl/pt;
829         tanth = TMath::Abs(tanth);
830         rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
831         ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
832
833         //   space-point resolutions
834         
835         sigmaRphi=SigmaY2(rpad,tanth,pt);
836         sigmaZ   =SigmaZ2(rpad,tanth   );
837         
838         //   cluster widths
839         
840         clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
841         clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
842         
843         // temporary protection
844         
845         if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
846         if(sigmaZ < 0.) sigmaZ=0.4e-3;
847         if(clRphi < 0.) clRphi=2.5e-3;
848         if(clZ < 0.) clZ=2.5e-5;
849         
850         //
851         
852         //
853         // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
854         // then the inaccuracy in a X-Y plane is only along Y (pad row)!
855         //
856         Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
857         Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
858         xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi));   // y
859           Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
860           fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
861           Float_t ymax=xprim*TMath::Tan(0.5*alpha);
862           if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim; 
863         xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
864           if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z(); 
865         xyz[2]=sigmaRphi;                                     // fSigmaY2
866         xyz[3]=sigmaZ;                                        // fSigmaZ2
867         xyz[4]=tpcHit->fQ;                                    // q
868
869         AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
870         if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);     
871
872         Int_t tracks[3]={tpcHit->Track(), -1, -1};
873         AliTPCcluster cluster(tracks,xyz);
874
875         clrow->InsertCluster(&cluster); nclusters++;
876
877         tpcHit = (AliTPChit*)NextHit();
878         
879
880       } // end of loop over hits
881
882     }   // end of loop over tracks
883
884     Int_t nrows=fTPCParam->GetNRow(isec);
885     for (Int_t irow=0; irow<nrows; irow++) {
886         AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
887         if (!clrow) continue;
888         carray.StoreRow(isec,irow);
889         carray.ClearRow(isec,irow);
890     }
891
892   } // end of loop over sectors  
893
894   cerr<<"Number of made clusters : "<<nclusters<<"                        \n";
895
896   carray.GetTree()->Write();
897
898   savedir->cd(); //switch back to the input file
899   
900 } // end of function
901
902 //_________________________________________________________________
903 void AliTPC::Hits2ExactClustersSector(Int_t isec)
904 {
905   //--------------------------------------------------------
906   //calculate exact cross point of track and given pad row
907   //resulting values are expressed in "digit" coordinata
908   //--------------------------------------------------------
909
910   //-----------------------------------------------------------------
911   // Origin: Marian Ivanov  GSI Darmstadt, m.ivanov@gsi.de
912   //-----------------------------------------------------------------
913   //
914   if (fClustersArray==0){    
915     return;
916   }
917   //
918   TParticle *particle; // pointer to a given particle
919   AliTPChit *tpcHit; // pointer to a sigle TPC hit
920   TClonesArray *particles; //pointer to the particle list
921   Int_t sector,nhits;
922   Int_t ipart;
923   const Int_t kcmaxhits=30000;
924   TVector * xxxx = new TVector(kcmaxhits*4);
925   TVector & xxx = *xxxx;
926   Int_t maxhits = kcmaxhits;
927   //construct array for each padrow
928   for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++) 
929     fClustersArray->CreateRow(isec,i);
930   
931   //---------------------------------------------------------------
932   //  Get the access to the tracks 
933   //---------------------------------------------------------------
934   
935   TTree *tH = gAlice->TreeH();
936   Stat_t ntracks = tH->GetEntries();
937   particles=gAlice->Particles();
938   Int_t npart = particles->GetEntriesFast();
939     
940   //------------------------------------------------------------
941   // Loop over tracks
942   //------------------------------------------------------------
943   
944   for(Int_t track=0;track<ntracks;track++){
945     ResetHits();
946     tH->GetEvent(track);
947     //
948     //  Get number of the TPC hits and a pointer
949     //  to the particles
950     //
951     nhits=fHits->GetEntriesFast();
952     //
953     // Loop over hits
954     //
955     Int_t currentIndex=0;
956     Int_t lastrow=-1;  //last writen row
957     for(Int_t hit=0;hit<nhits;hit++){
958       tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
959       if (tpcHit==0) continue;
960       sector=tpcHit->fSector; // sector number
961       if(sector != isec) continue; 
962       ipart=tpcHit->Track();
963       if (ipart<npart) particle=(TParticle*)particles->UncheckedAt(ipart);
964       
965       //find row number
966
967       Float_t  x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
968       Int_t    index[3]={1,isec,0};
969       Int_t    currentrow = fTPCParam->GetPadRow(x,index) ;     
970       if (currentrow<0) continue;
971       if (lastrow<0) lastrow=currentrow;
972       if (currentrow==lastrow){
973         if ( currentIndex>=maxhits){
974           maxhits+=kcmaxhits;
975           xxx.ResizeTo(4*maxhits);
976         }     
977         xxx(currentIndex*4)=x[0];
978         xxx(currentIndex*4+1)=x[1];
979         xxx(currentIndex*4+2)=x[2];     
980         xxx(currentIndex*4+3)=tpcHit->fQ;
981         currentIndex++; 
982       }
983       else 
984         if (currentIndex>2){
985           Float_t sumx=0;
986           Float_t sumx2=0;
987           Float_t sumx3=0;
988           Float_t sumx4=0;
989           Float_t sumy=0;
990           Float_t sumxy=0;
991           Float_t sumx2y=0;
992           Float_t sumz=0;
993           Float_t sumxz=0;
994           Float_t sumx2z=0;
995           Float_t sumq=0;
996           for (Int_t index=0;index<currentIndex;index++){
997             Float_t x,x2,x3,x4;
998             x=x2=x3=x4=xxx(index*4);
999             x2*=x;
1000             x3*=x2;
1001             x4*=x3;
1002             sumx+=x;
1003             sumx2+=x2;
1004             sumx3+=x3;
1005             sumx4+=x4;
1006             sumy+=xxx(index*4+1);
1007             sumxy+=xxx(index*4+1)*x;
1008             sumx2y+=xxx(index*4+1)*x2;
1009             sumz+=xxx(index*4+2);
1010             sumxz+=xxx(index*4+2)*x;
1011             sumx2z+=xxx(index*4+2)*x2;   
1012             sumq+=xxx(index*4+3);
1013           }
1014           Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1015           Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1016             sumx2*(sumx*sumx3-sumx2*sumx2);
1017           
1018           Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1019             sumx2*(sumxy*sumx3-sumx2y*sumx2);
1020           Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1021             sumx2*(sumxz*sumx3-sumx2z*sumx2);
1022           
1023           Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1024             sumx2*(sumx*sumx2y-sumx2*sumxy);
1025           Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1026             sumx2*(sumx*sumx2z-sumx2*sumxz);
1027           
1028           Float_t y=detay/det+centralPad;
1029           Float_t z=detaz/det;  
1030           Float_t by=detby/det; //y angle
1031           Float_t bz=detbz/det; //z angle
1032           sumy/=Float_t(currentIndex);
1033           sumz/=Float_t(currentIndex);
1034           AliComplexCluster cl;
1035           cl.fX=z;
1036           cl.fY=y;
1037           cl.fQ=sumq;
1038           cl.fSigmaX2=bz;
1039           cl.fSigmaY2=by;
1040           cl.fTracks[0]=ipart;
1041           
1042           AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1043           if (row!=0) row->InsertCluster(&cl);
1044           currentIndex=0;
1045           lastrow=currentrow;
1046         } //end of calculating cluster for given row
1047         
1048         
1049         
1050     } // end of loop over hits
1051   }   // end of loop over tracks 
1052   //write padrows to tree 
1053   for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1054     fClustersArray->StoreRow(isec,ii);    
1055     fClustersArray->ClearRow(isec,ii);        
1056   }
1057   xxxx->Delete();
1058  
1059 }
1060
1061 //__________________________________________________________________  
1062 void AliTPC::Hits2Digits()  
1063
1064  //----------------------------------------------------
1065  // Loop over all sectors
1066  //----------------------------------------------------
1067
1068   if(fTPCParam == 0){
1069     printf("AliTPCParam MUST be created firstly\n");
1070     return;
1071   } 
1072
1073  for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1074
1075 }
1076
1077
1078 //_____________________________________________________________________________
1079 void AliTPC::Hits2DigitsSector(Int_t isec)
1080 {
1081   //-------------------------------------------------------------------
1082   // TPC conversion from hits to digits.
1083   //------------------------------------------------------------------- 
1084
1085   //-----------------------------------------------------------------
1086   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1087   //-----------------------------------------------------------------
1088
1089   //-------------------------------------------------------
1090   //  Get the access to the track hits
1091   //-------------------------------------------------------
1092
1093
1094   TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1095   Stat_t ntracks = tH->GetEntries();
1096
1097   if( ntracks > 0){
1098
1099   //------------------------------------------- 
1100   //  Only if there are any tracks...
1101   //-------------------------------------------
1102
1103     TObjArray **row;
1104     
1105       printf("*** Processing sector number %d ***\n",isec);
1106
1107       Int_t nrows =fTPCParam->GetNRow(isec);
1108
1109       row= new TObjArray* [nrows];
1110     
1111       MakeSector(isec,nrows,tH,ntracks,row);
1112
1113       //--------------------------------------------------------
1114       //   Digitize this sector, row by row
1115       //   row[i] is the pointer to the TObjArray of TVectors,
1116       //   each one containing electrons accepted on this
1117       //   row, assigned into tracks
1118       //--------------------------------------------------------
1119
1120       Int_t i;
1121
1122       if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree();
1123
1124       for (i=0;i<nrows;i++){
1125
1126         AliDigits * dig = fDigitsArray->CreateRow(isec,i); 
1127
1128         DigitizeRow(i,isec,row);
1129
1130         fDigitsArray->StoreRow(isec,i);
1131
1132         Int_t ndig = dig->GetDigitSize(); 
1133  
1134         printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1135         
1136         fDigitsArray->ClearRow(isec,i);  
1137
1138    
1139        } // end of the sector digitization
1140
1141       for(i=0;i<nrows;i++){
1142         row[i]->Delete();  
1143         delete row[i];   
1144       }
1145       
1146        delete [] row; // delete the array of pointers to TObjArray-s
1147         
1148   } // ntracks >0
1149
1150 } // end of Hits2DigitsSector
1151
1152
1153 //_____________________________________________________________________________
1154 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1155 {
1156   //-----------------------------------------------------------
1157   // Single row digitization, coupling from the neighbouring
1158   // rows taken into account
1159   //-----------------------------------------------------------
1160
1161   //-----------------------------------------------------------------
1162   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1163   // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1164   //-----------------------------------------------------------------
1165  
1166
1167   Float_t zerosup = fTPCParam->GetZeroSup();
1168   Int_t nrows =fTPCParam->GetNRow(isec);
1169   fCurrentIndex[1]= isec;
1170   
1171
1172   Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1173   Int_t nofTbins = fTPCParam->GetMaxTBin();
1174   Int_t indexRange[4];
1175   //
1176   //  Integrated signal for this row
1177   //  and a single track signal
1178   //    
1179   TMatrix *m1   = new TMatrix(0,nofPads,0,nofTbins); // integrated
1180   TMatrix *m2   = new TMatrix(0,nofPads,0,nofTbins); // single
1181   //
1182   TMatrix &total  = *m1;
1183
1184   //  Array of pointers to the label-signal list
1185
1186   Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1187   Float_t  **pList = new Float_t* [nofDigits]; 
1188
1189   Int_t lp;
1190   Int_t i1;   
1191   for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1192   //
1193   //calculate signal 
1194   //
1195   Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1196   Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1197   for (Int_t row= row1;row<=row2;row++){
1198     Int_t nTracks= rows[row]->GetEntries();
1199     for (i1=0;i1<nTracks;i1++){
1200       fCurrentIndex[2]= row;
1201       fCurrentIndex[3]=irow;
1202       if (row==irow){
1203         m2->Zero();  // clear single track signal matrix
1204         Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange); 
1205         GetList(trackLabel,nofPads,m2,indexRange,pList);
1206       }
1207       else   GetSignal(rows[row],i1,0,m1,indexRange);
1208     }
1209   }
1210          
1211   Int_t tracks[3];
1212
1213   AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1214   for(Int_t ip=0;ip<nofPads;ip++){
1215     for(Int_t it=0;it<nofTbins;it++){
1216
1217       Float_t q = total(ip,it);
1218
1219       Int_t gi =it*nofPads+ip; // global index
1220
1221       q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac()); 
1222
1223       q = (Int_t)q;
1224
1225       if(q <=zerosup) continue; // do not fill zeros
1226       if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat();  // saturation
1227
1228       //
1229       //  "real" signal or electronic noise (list = -1)?
1230       //    
1231
1232       for(Int_t j1=0;j1<3;j1++){
1233         tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1234       }
1235
1236 //Begin_Html
1237 /*
1238   <A NAME="AliDigits"></A>
1239   using of AliDigits object
1240 */
1241 //End_Html
1242       dig->SetDigitFast((Short_t)q,it,ip);
1243       if (fDigitsArray->IsSimulated())
1244         {
1245          ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1246          ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1247          ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1248         }
1249      
1250     
1251     } // end of loop over time buckets
1252   }  // end of lop over pads 
1253
1254   //
1255   //  This row has been digitized, delete nonused stuff
1256   //
1257
1258   for(lp=0;lp<nofDigits;lp++){
1259     if(pList[lp]) delete [] pList[lp];
1260   }
1261   
1262   delete [] pList;
1263
1264   delete m1;
1265   delete m2;
1266   //  delete m3;
1267
1268 } // end of DigitizeRow
1269
1270 //_____________________________________________________________________________
1271
1272 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, TMatrix *m1, TMatrix *m2,
1273                           Int_t *indexRange)
1274 {
1275
1276   //---------------------------------------------------------------
1277   //  Calculates 2-D signal (pad,time) for a single track,
1278   //  returns a pointer to the signal matrix and the track label 
1279   //  No digitization is performed at this level!!!
1280   //---------------------------------------------------------------
1281
1282   //-----------------------------------------------------------------
1283   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1284   // Modified: Marian Ivanov 
1285   //-----------------------------------------------------------------
1286
1287   TVector *tv;
1288  
1289   tv = (TVector*)p1->At(ntr); // pointer to a track
1290   TVector &v = *tv;
1291   
1292   Float_t label = v(0);
1293   Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1294
1295   Int_t nElectrons = (tv->GetNrows()-1)/4;
1296   indexRange[0]=9999; // min pad
1297   indexRange[1]=-1; // max pad
1298   indexRange[2]=9999; //min time
1299   indexRange[3]=-1; // max time
1300
1301   //  Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1302
1303   TMatrix &signal = *m1;
1304   TMatrix &total = *m2;
1305   //
1306   //  Loop over all electrons
1307   //
1308   for(Int_t nel=0; nel<nElectrons; nel++){
1309     Int_t idx=nel*4;
1310     Float_t aval =  v(idx+4);
1311     Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac(); 
1312     Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1313     Int_t n = fTPCParam->CalcResponse(xyz,fCurrentIndex,fCurrentIndex[3]);
1314     
1315     if (n>0) for (Int_t i =0; i<n; i++){
1316        Int_t *index = fTPCParam->GetResBin(i);        
1317        Int_t pad=index[1]+centralPad;  //in digit coordinates central pad has coordinate 0
1318        if ( ( pad<(fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]))) && (pad>0)) {
1319          Int_t time=index[2];    
1320          Float_t weight = fTPCParam->GetResWeight(i); //we normalise response to ADC channel
1321          weight *= eltoadcfac;
1322          
1323          if (m1!=0) signal(pad,time)+=weight; 
1324          total(pad,time)+=weight;
1325          indexRange[0]=TMath::Min(indexRange[0],pad);
1326          indexRange[1]=TMath::Max(indexRange[1],pad);
1327          indexRange[2]=TMath::Min(indexRange[2],time);
1328          indexRange[3]=TMath::Max(indexRange[3],time); 
1329        }         
1330     }
1331   } // end of loop over electrons
1332   
1333   return label; // returns track label when finished
1334 }
1335
1336 //_____________________________________________________________________________
1337 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *indexRange,
1338                      Float_t **pList)
1339 {
1340   //----------------------------------------------------------------------
1341   //  Updates the list of tracks contributing to digits for a given row
1342   //----------------------------------------------------------------------
1343
1344   //-----------------------------------------------------------------
1345   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1346   //-----------------------------------------------------------------
1347
1348   TMatrix &signal = *m;
1349
1350   // lop over nonzero digits
1351
1352   for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1353     for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1354
1355
1356         // accept only the contribution larger than 500 electrons (1/2 s_noise)
1357
1358         if(signal(ip,it)<0.5) continue; 
1359
1360
1361         Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1362         
1363         if(!pList[globalIndex]){
1364         
1365           // 
1366           // Create new list (6 elements - 3 signals and 3 labels),
1367           //
1368
1369           pList[globalIndex] = new Float_t [6];
1370
1371           // set list to -1 
1372
1373           *pList[globalIndex] = -1.;
1374           *(pList[globalIndex]+1) = -1.;
1375           *(pList[globalIndex]+2) = -1.;
1376           *(pList[globalIndex]+3) = -1.;
1377           *(pList[globalIndex]+4) = -1.;
1378           *(pList[globalIndex]+5) = -1.;
1379
1380
1381           *pList[globalIndex] = label;
1382           *(pList[globalIndex]+3) = signal(ip,it);
1383         }
1384         else{
1385
1386           // check the signal magnitude
1387
1388           Float_t highest = *(pList[globalIndex]+3);
1389           Float_t middle = *(pList[globalIndex]+4);
1390           Float_t lowest = *(pList[globalIndex]+5);
1391
1392           //
1393           //  compare the new signal with already existing list
1394           //
1395
1396           if(signal(ip,it)<lowest) continue; // neglect this track
1397
1398           //
1399
1400           if (signal(ip,it)>highest){
1401             *(pList[globalIndex]+5) = middle;
1402             *(pList[globalIndex]+4) = highest;
1403             *(pList[globalIndex]+3) = signal(ip,it);
1404
1405             *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1406             *(pList[globalIndex]+1) = *pList[globalIndex];
1407             *pList[globalIndex] = label;
1408           }
1409           else if (signal(ip,it)>middle){
1410             *(pList[globalIndex]+5) = middle;
1411             *(pList[globalIndex]+4) = signal(ip,it);
1412
1413             *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1414             *(pList[globalIndex]+1) = label;
1415           }
1416           else{
1417             *(pList[globalIndex]+5) = signal(ip,it);
1418             *(pList[globalIndex]+2) = label;
1419           }
1420         }
1421
1422     } // end of loop over pads
1423   } // end of loop over time bins
1424
1425
1426
1427 }//end of GetList
1428 //___________________________________________________________________
1429 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1430                         Stat_t ntracks,TObjArray **row)
1431 {
1432
1433   //-----------------------------------------------------------------
1434   // Prepares the sector digitization, creates the vectors of
1435   // tracks for each row of this sector. The track vector
1436   // contains the track label and the position of electrons.
1437   //-----------------------------------------------------------------
1438
1439   //-----------------------------------------------------------------
1440   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1441   //-----------------------------------------------------------------
1442
1443   Float_t gasgain = fTPCParam->GetGasGain();
1444   Int_t i;
1445   Float_t xyz[4]; 
1446
1447   AliTPChit *tpcHit; // pointer to a sigle TPC hit    
1448   //MI change
1449   TBranch * branch=0;
1450   if (fHitType&2) branch = TH->GetBranch("TPC2");
1451   else branch = TH->GetBranch("TPC");
1452
1453  
1454   //----------------------------------------------
1455   // Create TObjArray-s, one for each row,
1456   // each TObjArray will store the TVectors
1457   // of electrons, one TVector per each track.
1458   //---------------------------------------------- 
1459     
1460   Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
1461   TVector **tracks = new TVector* [nrows]; //pointers to the track vectors
1462   for(i=0; i<nrows; i++){
1463     row[i] = new TObjArray;
1464     nofElectrons[i]=0;
1465     tracks[i]=0;
1466   }
1467
1468  
1469
1470   //--------------------------------------------------------------------
1471   //  Loop over tracks, the "track" contains the full history
1472   //--------------------------------------------------------------------
1473
1474   Int_t previousTrack,currentTrack;
1475   previousTrack = -1; // nothing to store so far!
1476
1477   for(Int_t track=0;track<ntracks;track++){
1478     Bool_t isInSector=kTRUE;
1479     ResetHits();
1480
1481     if (fHitType&2) {
1482       isInSector=kFALSE;
1483       TBranch * br = TH->GetBranch("fTrackHitsInfo");
1484       br->GetEvent(track);
1485       AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
1486       for (UInt_t j=0;j<ar->GetSize();j++){
1487         if (  ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
1488       }
1489     }
1490     if (!isInSector) continue;
1491     //MI change
1492     branch->GetEntry(track); // get next track
1493
1494     //M.I. changes
1495
1496     tpcHit = (AliTPChit*)FirstHit(-1);
1497
1498     //--------------------------------------------------------------
1499     //  Loop over hits
1500     //--------------------------------------------------------------
1501
1502
1503     while(tpcHit){
1504       
1505       Int_t sector=tpcHit->fSector; // sector number
1506       //      if(sector != isec) continue; 
1507       if(sector != isec){
1508         tpcHit = (AliTPChit*) NextHit();
1509         continue; 
1510       }
1511
1512         currentTrack = tpcHit->Track(); // track number
1513
1514
1515         if(currentTrack != previousTrack){
1516                           
1517            // store already filled fTrack
1518               
1519            for(i=0;i<nrows;i++){
1520              if(previousTrack != -1){
1521                if(nofElectrons[i]>0){
1522                  TVector &v = *tracks[i];
1523                  v(0) = previousTrack;
1524                  tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1525                  row[i]->Add(tracks[i]);                     
1526                }
1527                else{
1528                  delete tracks[i]; // delete empty TVector
1529                  tracks[i]=0;
1530                }
1531              }
1532
1533              nofElectrons[i]=0;
1534              tracks[i] = new TVector(481); // TVectors for the next fTrack
1535
1536            } // end of loop over rows
1537                
1538            previousTrack=currentTrack; // update track label 
1539         }
1540            
1541         Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1542
1543        //---------------------------------------------------
1544        //  Calculate the electron attachment probability
1545        //---------------------------------------------------
1546
1547
1548         Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
1549                                                         /fTPCParam->GetDriftV(); 
1550         // in microseconds!     
1551         Float_t attProb = fTPCParam->GetAttCoef()*
1552           fTPCParam->GetOxyCont()*time; //  fraction! 
1553    
1554         //-----------------------------------------------
1555         //  Loop over electrons
1556         //-----------------------------------------------
1557         Int_t index[3];
1558         index[1]=isec;
1559         for(Int_t nel=0;nel<qI;nel++){
1560           // skip if electron lost due to the attachment
1561           if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1562           xyz[0]=tpcHit->X();
1563           xyz[1]=tpcHit->Y();
1564           xyz[2]=tpcHit->Z();     
1565           xyz[3]= (Float_t) (-gasgain*TMath::Log(gRandom->Rndm()));
1566           index[0]=1;
1567           
1568           TransportElectron(xyz,index); //MI change -august       
1569           Int_t rowNumber;
1570           fTPCParam->GetPadRow(xyz,index); //MI change august
1571           rowNumber = index[2];
1572           //transform position to local digit coordinates
1573           //relative to nearest pad row 
1574           if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;       
1575           nofElectrons[rowNumber]++;      
1576           //----------------------------------
1577           // Expand vector if necessary
1578           //----------------------------------
1579           if(nofElectrons[rowNumber]>120){
1580             Int_t range = tracks[rowNumber]->GetNrows();
1581             if((nofElectrons[rowNumber])>(range-1)/4){
1582         
1583               tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
1584             }
1585           }
1586           
1587           TVector &v = *tracks[rowNumber];
1588           Int_t idx = 4*nofElectrons[rowNumber]-3;
1589
1590           v(idx)=  xyz[0];   // X - pad row coordinate
1591           v(idx+1)=xyz[1];   // Y - pad coordinate (along the pad-row)
1592           v(idx+2)=xyz[2];   // Z - time bin coordinate
1593           v(idx+3)=xyz[3];   // avalanche size  
1594         } // end of loop over electrons
1595
1596         tpcHit = (AliTPChit*)NextHit();
1597         
1598       } // end of loop over hits
1599     } // end of loop over tracks
1600
1601     //
1602     //   store remaining track (the last one) if not empty
1603     //
1604
1605      for(i=0;i<nrows;i++){
1606        if(nofElectrons[i]>0){
1607           TVector &v = *tracks[i];
1608           v(0) = previousTrack;
1609           tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1610           row[i]->Add(tracks[i]);  
1611         }
1612         else{
1613           delete tracks[i];
1614           tracks[i]=0;
1615         }  
1616       }  
1617
1618           delete [] tracks;
1619           delete [] nofElectrons;
1620  
1621
1622 } // end of MakeSector
1623
1624
1625 //_____________________________________________________________________________
1626 void AliTPC::Init()
1627 {
1628   //
1629   // Initialise TPC detector after definition of geometry
1630   //
1631   Int_t i;
1632   //
1633   printf("\n");
1634   for(i=0;i<35;i++) printf("*");
1635   printf(" TPC_INIT ");
1636   for(i=0;i<35;i++) printf("*");
1637   printf("\n");
1638   //
1639   for(i=0;i<80;i++) printf("*");
1640   printf("\n");
1641 }
1642
1643 //_____________________________________________________________________________
1644 void AliTPC::MakeBranch(Option_t* option)
1645 {
1646   //
1647   // Create Tree branches for the TPC.
1648   //
1649   Int_t buffersize = 4000;
1650   char branchname[10];
1651   sprintf(branchname,"%s",GetName());
1652
1653   AliDetector::MakeBranch(option);
1654
1655   char *d = strstr(option,"D");
1656
1657   if (fDigits   && gAlice->TreeD() && d) {
1658     gAlice->TreeD()->Branch(branchname,&fDigits, buffersize);
1659     printf("Making Branch %s for digits\n",branchname);
1660   }     
1661
1662   if (fHitType&2) MakeBranch2(option); // MI change 14.09.2000
1663 }
1664  
1665 //_____________________________________________________________________________
1666 void AliTPC::ResetDigits()
1667 {
1668   //
1669   // Reset number of digits and the digits array for this detector
1670   //
1671   fNdigits   = 0;
1672   if (fDigits)   fDigits->Clear();
1673 }
1674
1675 //_____________________________________________________________________________
1676 void AliTPC::SetSecAL(Int_t sec)
1677 {
1678   //---------------------------------------------------
1679   // Activate/deactivate selection for lower sectors
1680   //---------------------------------------------------
1681
1682   //-----------------------------------------------------------------
1683   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1684   //-----------------------------------------------------------------
1685
1686   fSecAL = sec;
1687 }
1688
1689 //_____________________________________________________________________________
1690 void AliTPC::SetSecAU(Int_t sec)
1691 {
1692   //----------------------------------------------------
1693   // Activate/deactivate selection for upper sectors
1694   //---------------------------------------------------
1695
1696   //-----------------------------------------------------------------
1697   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1698   //-----------------------------------------------------------------
1699
1700   fSecAU = sec;
1701 }
1702
1703 //_____________________________________________________________________________
1704 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
1705 {
1706   //----------------------------------------
1707   // Select active lower sectors
1708   //----------------------------------------
1709
1710   //-----------------------------------------------------------------
1711   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1712   //-----------------------------------------------------------------
1713
1714   fSecLows[0] = s1;
1715   fSecLows[1] = s2;
1716   fSecLows[2] = s3;
1717   fSecLows[3] = s4;
1718   fSecLows[4] = s5;
1719   fSecLows[5] = s6;
1720 }
1721
1722 //_____________________________________________________________________________
1723 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
1724                        Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10, 
1725                        Int_t s11 , Int_t s12)
1726 {
1727   //--------------------------------
1728   // Select active upper sectors
1729   //--------------------------------
1730
1731   //-----------------------------------------------------------------
1732   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1733   //-----------------------------------------------------------------
1734
1735   fSecUps[0] = s1;
1736   fSecUps[1] = s2;
1737   fSecUps[2] = s3;
1738   fSecUps[3] = s4;
1739   fSecUps[4] = s5;
1740   fSecUps[5] = s6;
1741   fSecUps[6] = s7;
1742   fSecUps[7] = s8;
1743   fSecUps[8] = s9;
1744   fSecUps[9] = s10;
1745   fSecUps[10] = s11;
1746   fSecUps[11] = s12;
1747 }
1748
1749 //_____________________________________________________________________________
1750 void AliTPC::SetSens(Int_t sens)
1751 {
1752
1753   //-------------------------------------------------------------
1754   // Activates/deactivates the sensitive strips at the center of
1755   // the pad row -- this is for the space-point resolution calculations
1756   //-------------------------------------------------------------
1757
1758   //-----------------------------------------------------------------
1759   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1760   //-----------------------------------------------------------------
1761
1762   fSens = sens;
1763 }
1764
1765  
1766 void AliTPC::SetSide(Float_t side=0.)
1767 {
1768   // choice of the TPC side
1769
1770   fSide = side;
1771  
1772 }
1773 //____________________________________________________________________________
1774 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
1775                            Float_t p2,Float_t p3)
1776 {
1777
1778   // gax mixture definition
1779
1780  fNoComp = nc;
1781  
1782  fMixtComp[0]=c1;
1783  fMixtComp[1]=c2;
1784  fMixtComp[2]=c3;
1785
1786  fMixtProp[0]=p1;
1787  fMixtProp[1]=p2;
1788  fMixtProp[2]=p3; 
1789  
1790  
1791 }
1792 //_____________________________________________________________________________
1793
1794 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
1795 {
1796   //
1797   // electron transport taking into account:
1798   // 1. diffusion, 
1799   // 2.ExB at the wires
1800   // 3. nonisochronity
1801   //
1802   // xyz and index must be already transformed to system 1
1803   //
1804
1805   fTPCParam->Transform1to2(xyz,index);
1806   
1807   //add diffusion
1808   Float_t driftl=xyz[2];
1809   if(driftl<0.01) driftl=0.01;
1810   driftl=TMath::Sqrt(driftl);
1811   Float_t sigT = driftl*(fTPCParam->GetDiffT());
1812   Float_t sigL = driftl*(fTPCParam->GetDiffL());
1813   xyz[0]=gRandom->Gaus(xyz[0],sigT);
1814   xyz[1]=gRandom->Gaus(xyz[1],sigT);
1815   xyz[2]=gRandom->Gaus(xyz[2],sigL);
1816
1817   // ExB
1818   
1819   if (fTPCParam->GetMWPCReadout()==kTRUE){
1820     Float_t x1=xyz[0];
1821     fTPCParam->Transform2to2NearestWire(xyz,index);
1822     Float_t dx=xyz[0]-x1;
1823     xyz[1]+=dx*(fTPCParam->GetOmegaTau());
1824   }
1825   //add nonisochronity (not implemented yet)
1826   
1827 }
1828 //_____________________________________________________________________________
1829 void AliTPC::Streamer(TBuffer &R__b)
1830 {
1831   //
1832   // Stream an object of class AliTPC.
1833   //
1834    if (R__b.IsReading()) {
1835       Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1836       AliDetector::Streamer(R__b);
1837       R__b >> fTPCParam;
1838       if (R__v < 2) return;
1839       R__b >> fNsectors;
1840       R__b >> fHitType;
1841
1842    } else {
1843       R__b.WriteVersion(AliTPC::IsA());
1844       AliDetector::Streamer(R__b);
1845       R__b << fTPCParam;
1846       R__b << fNsectors;
1847       R__b << fHitType; 
1848    }
1849 }
1850   
1851 ClassImp(AliTPCdigit)
1852  
1853 //_____________________________________________________________________________
1854 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
1855   AliDigit(tracks)
1856 {
1857   //
1858   // Creates a TPC digit object
1859   //
1860   fSector     = digits[0];
1861   fPadRow     = digits[1];
1862   fPad        = digits[2];
1863   fTime       = digits[3];
1864   fSignal     = digits[4];
1865 }
1866
1867  
1868 ClassImp(AliTPChit)
1869  
1870 //_____________________________________________________________________________
1871 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
1872 AliHit(shunt,track)
1873 {
1874   //
1875   // Creates a TPC hit object
1876   //
1877   fSector     = vol[0];
1878   fPadRow     = vol[1];
1879   fX          = hits[0];
1880   fY          = hits[1];
1881   fZ          = hits[2];
1882   fQ          = hits[3];
1883 }
1884  
1885
1886 //________________________________________________________________________
1887 // Additional code because of the AliTPCTrackHits
1888
1889 void AliTPC::MakeBranch2(Option_t *option)
1890 {
1891   //
1892   // Create a new branch in the current Root Tree
1893   // The branch of fHits is automatically split
1894   // MI change 14.09.2000
1895   if (fHitType&2==0) return;
1896   char branchname[10];
1897   sprintf(branchname,"%s2",GetName());  
1898   //
1899   // Get the pointer to the header
1900   char *cH = strstr(option,"H");
1901   //
1902   if (fTrackHits   && gAlice->TreeH() && cH) {    
1903     //    gAlice->TreeH()->Branch(branchname,&fTrackHits, fBufferSize);
1904     AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHits, 
1905                                                    gAlice->TreeH(),fBufferSize,1);
1906     gAlice->TreeH()->GetListOfBranches()->Add(branch);
1907     printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
1908   }     
1909 }
1910
1911 void AliTPC::SetTreeAddress()
1912 {
1913   if (fHitType&1) AliDetector::SetTreeAddress();
1914   if (fHitType&2) SetTreeAddress2();
1915 }
1916
1917 void AliTPC::SetTreeAddress2()
1918 {
1919   //
1920   // Set branch address for the TrackHits Tree
1921   // 
1922   TBranch *branch;
1923   char branchname[20];
1924   sprintf(branchname,"%s2",GetName());
1925   //
1926   // Branch address for hit tree
1927   TTree *treeH = gAlice->TreeH();
1928   if (treeH) {
1929     branch = treeH->GetBranch(branchname);
1930     if (branch) branch->SetAddress(&fTrackHits);
1931   }
1932 }
1933
1934 void AliTPC::FinishPrimary()
1935 {
1936   if (fTrackHits) fTrackHits->FlushHitStack();  
1937 }
1938
1939
1940 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
1941
1942   //
1943   // add hit to the list  
1944   TClonesArray &particles = *(gAlice->Particles());
1945   Int_t rtrack;
1946   if (fIshunt) {
1947     int primary = gAlice->GetPrimary(track);
1948     ((TParticle *)particles[primary])->SetBit(kKeepBit);
1949     rtrack=primary;
1950   } else {
1951     rtrack=track;
1952     gAlice->FlagTrack(track);
1953   }  
1954   //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
1955   //if (hit->fTrack!=rtrack)
1956   //  cout<<"bad track number\n";
1957   if (fTrackHits) 
1958     fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
1959                              hits[1],hits[2],(Int_t)hits[3]);
1960 }
1961
1962 void AliTPC::ResetHits()
1963 {
1964   if (fHitType&1) AliDetector::ResetHits();
1965   if (fHitType&2) ResetHits2();
1966 }
1967
1968 void AliTPC::ResetHits2()
1969 {
1970   //
1971   //reset hits
1972   if (fTrackHits) fTrackHits->Clear();
1973 }   
1974
1975 AliHit* AliTPC::FirstHit(Int_t track)
1976 {
1977   if (fHitType&2) return FirstHit2(track);
1978   return AliDetector::FirstHit(track);
1979 }
1980 AliHit* AliTPC::NextHit()
1981 {
1982   if (fHitType&2) return NextHit2();
1983   return AliDetector::NextHit();
1984 }
1985
1986 AliHit* AliTPC::FirstHit2(Int_t track)
1987 {
1988   //
1989   // Initialise the hit iterator
1990   // Return the address of the first hit for track
1991   // If track>=0 the track is read from disk
1992   // while if track<0 the first hit of the current
1993   // track is returned
1994   // 
1995   if(track>=0) {
1996     gAlice->ResetHits();
1997     gAlice->TreeH()->GetEvent(track);
1998   }
1999   //
2000   if (fTrackHits) {
2001     fTrackHits->First();
2002     return fTrackHits->GetHit();
2003   }
2004   else return 0;
2005 }
2006
2007 AliHit* AliTPC::NextHit2()
2008 {
2009   //
2010   //Return the next hit for the current track
2011
2012   if (fTrackHits) {
2013     fTrackHits->Next();
2014     return fTrackHits->GetHit();
2015   }
2016   else 
2017     return 0;
2018 }
2019
2020 void AliTPC::LoadPoints(Int_t)
2021 {
2022   //
2023   Int_t a = 0;
2024   LoadPoints2(a);
2025   // LoadPoints3(a);
2026
2027 }
2028
2029
2030 void AliTPC::RemapTrackHitIDs(Int_t *map)
2031 {
2032   if (!fTrackHits) return;
2033   AliObjectArray * arr = fTrackHits->fTrackHitsInfo;
2034   for (UInt_t i=0;i<arr->GetSize();i++){
2035     AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2036     info->fTrackID = map[info->fTrackID];
2037   }
2038   
2039 }
2040
2041
2042 //_____________________________________________________________________________
2043 void AliTPC::LoadPoints2(Int_t)
2044 {
2045   //
2046   // Store x, y, z of all hits in memory
2047   //
2048   if (fTrackHits == 0) return;
2049   //
2050   Int_t nhits = fTrackHits->GetEntriesFast();
2051   if (nhits == 0) return;
2052   Int_t tracks = gAlice->GetNtrack();
2053   if (fPoints == 0) fPoints = new TObjArray(tracks);
2054   AliHit *ahit;
2055   //
2056   Int_t *ntrk=new Int_t[tracks];
2057   Int_t *limi=new Int_t[tracks];
2058   Float_t **coor=new Float_t*[tracks];
2059   for(Int_t i=0;i<tracks;i++) {
2060     ntrk[i]=0;
2061     coor[i]=0;
2062     limi[i]=0;
2063   }
2064   //
2065   AliPoints *points = 0;
2066   Float_t *fp=0;
2067   Int_t trk;
2068   Int_t chunk=nhits/4+1;
2069   //
2070   // Loop over all the hits and store their position
2071   //
2072   ahit = FirstHit2(-1);
2073   //for (Int_t hit=0;hit<nhits;hit++) {
2074   while (ahit){
2075     //    ahit = (AliHit*)fHits->UncheckedAt(hit);
2076     trk=ahit->GetTrack();
2077     if(ntrk[trk]==limi[trk]) {
2078       //
2079       // Initialise a new track
2080       fp=new Float_t[3*(limi[trk]+chunk)];
2081       if(coor[trk]) {
2082         memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2083         delete [] coor[trk];
2084       }
2085       limi[trk]+=chunk;
2086       coor[trk] = fp;
2087     } else {
2088       fp = coor[trk];
2089     }
2090     fp[3*ntrk[trk]  ] = ahit->X();
2091     fp[3*ntrk[trk]+1] = ahit->Y();
2092     fp[3*ntrk[trk]+2] = ahit->Z();
2093     ntrk[trk]++;
2094     ahit = NextHit2();
2095   }
2096   //
2097   for(trk=0; trk<tracks; ++trk) {
2098     if(ntrk[trk]) {
2099       points = new AliPoints();
2100       points->SetMarkerColor(GetMarkerColor());
2101       points->SetMarkerSize(GetMarkerSize());
2102       points->SetDetector(this);
2103       points->SetParticle(trk);
2104       points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2105       fPoints->AddAt(points,trk);
2106       delete [] coor[trk];
2107       coor[trk]=0;
2108     }
2109   }
2110   delete [] coor;
2111   delete [] ntrk;
2112   delete [] limi;
2113 }
2114
2115
2116 //_____________________________________________________________________________
2117 void AliTPC::LoadPoints3(Int_t)
2118 {
2119   //
2120   // Store x, y, z of all hits in memory
2121   // - only intersection point with pad row
2122   if (fTrackHits == 0) return;
2123   //
2124   Int_t nhits = fTrackHits->GetEntriesFast();
2125   if (nhits == 0) return;
2126   Int_t tracks = gAlice->GetNtrack();
2127   if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2128   fPoints->Expand(2*tracks);
2129   AliHit *ahit;
2130   //
2131   Int_t *ntrk=new Int_t[tracks];
2132   Int_t *limi=new Int_t[tracks];
2133   Float_t **coor=new Float_t*[tracks];
2134   for(Int_t i=0;i<tracks;i++) {
2135     ntrk[i]=0;
2136     coor[i]=0;
2137     limi[i]=0;
2138   }
2139   //
2140   AliPoints *points = 0;
2141   Float_t *fp=0;
2142   Int_t trk;
2143   Int_t chunk=nhits/4+1;
2144   //
2145   // Loop over all the hits and store their position
2146   //
2147   ahit = FirstHit2(-1);
2148   //for (Int_t hit=0;hit<nhits;hit++) {
2149
2150   Int_t lastrow = -1;
2151   while (ahit){
2152     //    ahit = (AliHit*)fHits->UncheckedAt(hit);
2153     trk=ahit->GetTrack(); 
2154     Float_t  x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2155     Int_t    index[3]={1,((AliTPChit*)ahit)->fSector,0};
2156     Int_t    currentrow = fTPCParam->GetPadRow(x,index) ;
2157     if (currentrow!=lastrow){
2158       lastrow = currentrow;
2159       //later calculate intersection point           
2160       if(ntrk[trk]==limi[trk]) {
2161         //
2162         // Initialise a new track
2163         fp=new Float_t[3*(limi[trk]+chunk)];
2164         if(coor[trk]) {
2165           memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2166           delete [] coor[trk];
2167         }
2168         limi[trk]+=chunk;
2169         coor[trk] = fp;
2170       } else {
2171         fp = coor[trk];
2172       }
2173       fp[3*ntrk[trk]  ] = ahit->X();
2174       fp[3*ntrk[trk]+1] = ahit->Y();
2175       fp[3*ntrk[trk]+2] = ahit->Z();
2176       ntrk[trk]++;
2177     }
2178     ahit = NextHit2();
2179   }
2180   
2181   //
2182   for(trk=0; trk<tracks; ++trk) {
2183     if(ntrk[trk]) {
2184       points = new AliPoints();
2185       points->SetMarkerColor(GetMarkerColor()+1);
2186       points->SetMarkerStyle(5);
2187       points->SetMarkerSize(0.2);
2188       points->SetDetector(this);
2189       points->SetParticle(trk);
2190       //      points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2191       points->SetPolyMarker(ntrk[trk],coor[trk],30);
2192       fPoints->AddAt(points,tracks+trk);
2193       delete [] coor[trk];
2194       coor[trk]=0;
2195     }
2196   }
2197   delete [] coor;
2198   delete [] ntrk;
2199   delete [] limi;
2200 }
2201
2202
2203
2204 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2205 {
2206
2207   //
2208   //fill clones array with intersection of current point with the
2209   //middle of the row
2210   Int_t sector;
2211   Int_t ipart;
2212   
2213   const Int_t kcmaxhits=30000;
2214   TVector * xxxx = new TVector(kcmaxhits*4);
2215   TVector & xxx = *xxxx;
2216   Int_t maxhits = kcmaxhits;
2217       
2218   //
2219   AliTPChit * tpcHit=0;
2220   tpcHit = (AliTPChit*)FirstHit2(-1);
2221   Int_t currentIndex=0;
2222   Int_t lastrow=-1;  //last writen row
2223
2224   while (tpcHit){
2225     if (tpcHit==0) continue;
2226     sector=tpcHit->fSector; // sector number
2227     ipart=tpcHit->Track();
2228     
2229     //find row number
2230     
2231     Float_t  x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2232     Int_t    index[3]={1,sector,0};
2233     Int_t    currentrow = fTPCParam->GetPadRow(x,index) ;       
2234     if (currentrow<0) continue;
2235     if (lastrow<0) lastrow=currentrow;
2236     if (currentrow==lastrow){
2237       if ( currentIndex>=maxhits){
2238         maxhits+=kcmaxhits;
2239         xxx.ResizeTo(4*maxhits);
2240       }     
2241       xxx(currentIndex*4)=x[0];
2242       xxx(currentIndex*4+1)=x[1];
2243       xxx(currentIndex*4+2)=x[2];       
2244       xxx(currentIndex*4+3)=tpcHit->fQ;
2245       currentIndex++;   
2246     }
2247     else 
2248       if (currentIndex>2){
2249         Float_t sumx=0;
2250         Float_t sumx2=0;
2251         Float_t sumx3=0;
2252         Float_t sumx4=0;
2253         Float_t sumy=0;
2254         Float_t sumxy=0;
2255         Float_t sumx2y=0;
2256         Float_t sumz=0;
2257         Float_t sumxz=0;
2258         Float_t sumx2z=0;
2259         Float_t sumq=0;
2260         for (Int_t index=0;index<currentIndex;index++){
2261           Float_t x,x2,x3,x4;
2262           x=x2=x3=x4=xxx(index*4);
2263           x2*=x;
2264           x3*=x2;
2265           x4*=x3;
2266           sumx+=x;
2267           sumx2+=x2;
2268           sumx3+=x3;
2269           sumx4+=x4;
2270           sumy+=xxx(index*4+1);
2271           sumxy+=xxx(index*4+1)*x;
2272           sumx2y+=xxx(index*4+1)*x2;
2273           sumz+=xxx(index*4+2);
2274           sumxz+=xxx(index*4+2)*x;
2275           sumx2z+=xxx(index*4+2)*x2;     
2276           sumq+=xxx(index*4+3);
2277         }
2278         Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2279         Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2280           sumx2*(sumx*sumx3-sumx2*sumx2);
2281         
2282         Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2283           sumx2*(sumxy*sumx3-sumx2y*sumx2);
2284         Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2285           sumx2*(sumxz*sumx3-sumx2z*sumx2);
2286         
2287         Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2288           sumx2*(sumx*sumx2y-sumx2*sumxy);
2289         Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2290           sumx2*(sumx*sumx2z-sumx2*sumxz);
2291         
2292         Float_t y=detay/det+centralPad;
2293         Float_t z=detaz/det;    
2294         Float_t by=detby/det; //y angle
2295         Float_t bz=detbz/det; //z angle
2296         sumy/=Float_t(currentIndex);
2297         sumz/=Float_t(currentIndex);
2298         
2299         AliComplexCluster cl;
2300         cl.fX=z;
2301         cl.fY=y;
2302         cl.fQ=sumq;
2303         cl.fSigmaX2=bz;
2304         cl.fSigmaY2=by;
2305         cl.fTracks[0]=ipart;
2306         
2307         AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2308         if (row!=0) row->InsertCluster(&cl);
2309         currentIndex=0;
2310         lastrow=currentrow;
2311       } //end of calculating cluster for given row
2312                 
2313   } // end of loop over hits
2314   xxxx->Delete();
2315  
2316
2317
2318 }