]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPC.cxx
Transition to NewIO
[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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  Time Projection Chamber                                                  //
21 //  This class contains the basic functions for the Time Projection Chamber  //
22 //  detector. Functions specific to one particular geometry are              //
23 //  contained in the derived classes                                         //
24 //                                                                           //
25 //Begin_Html
26 /*
27 <img src="picts/AliTPCClass.gif">
28 */
29 //End_Html
30 //                                                                           //
31 //                                                                          //
32 ///////////////////////////////////////////////////////////////////////////////
33
34 //
35
36 #include <Riostream.h>
37 #include <stdlib.h>
38
39 #include <TFile.h>  
40 #include <TGeometry.h>
41 #include <TInterpreter.h>
42 #include <TMath.h>
43 #include <TMatrix.h>
44 #include <TNode.h>
45 #include <TObjectTable.h>
46 #include <TParticle.h>
47 #include <TROOT.h>
48 #include <TRandom.h>
49 #include <TSystem.h>     
50 #include <TTUBS.h>
51 #include <TTree.h>
52 #include <TVector.h>
53 #include <TVirtualMC.h>
54
55 #include "AliArrayBranch.h"
56 #include "AliClusters.h"
57 #include "AliComplexCluster.h"
58 #include "AliDigits.h"
59 #include "AliMagF.h"
60 #include "AliPoints.h"
61 #include "AliRun.h"
62 #include "AliRunLoader.h"
63 #include "AliSimDigits.h"
64 #include "AliTPC.h"
65 #include "AliTPC.h"
66 #include "AliTPCClustersArray.h"
67 #include "AliTPCClustersRow.h"
68 #include "AliTPCDigitsArray.h"
69 #include "AliTPCLoader.h"
70 #include "AliTPCPRF2D.h"
71 #include "AliTPCParamSR.h"
72 #include "AliTPCRF1D.h"
73 #include "AliTPCTrackHits.h"
74 #include "AliTPCTrackHitsV2.h"
75 #include "AliTPCcluster.h"
76 #include "AliTPCclusterer.h"
77 #include "AliTPCtracker.h"
78 #include "AliTrackReference.h"
79
80
81 ClassImp(AliTPC) 
82
83 //_____________________________________________________________________________
84 // helper class for fast matrix and vector manipulation - no range checking
85 // origin - Marian Ivanov
86
87 class AliTPCFastMatrix : public TMatrix {
88 public :
89   AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb);
90   inline Float_t & UncheckedAt(Int_t rown, Int_t coln) const  {return  (fIndex[coln])[rown];} //fast acces   
91   inline Float_t   UncheckedAtFast(Int_t rown, Int_t coln) const  {return  (fIndex[coln])[rown];} //fast acces   
92 };
93
94 AliTPCFastMatrix::AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb):
95   TMatrix(row_lwb, row_upb,col_lwb,col_upb)
96    {
97    };
98 //_____________________________________________________________________________
99 class AliTPCFastVector : public TVector {
100 public :
101   AliTPCFastVector(Int_t size);
102   virtual ~AliTPCFastVector(){;}
103   inline Float_t & UncheckedAt(Int_t index) const  {return  fElements[index];} //fast acces  
104   
105 };
106
107 AliTPCFastVector::AliTPCFastVector(Int_t size):TVector(size){
108 };
109
110 //_____________________________________________________________________________
111 AliTPC::AliTPC()
112 {
113   //
114   // Default constructor
115   //
116   fIshunt   = 0;
117   fHits     = 0;
118   fDigits   = 0;
119   fNsectors = 0;
120   fDigitsArray = 0;
121   fClustersArray = 0;
122   fDefaults = 0;
123   fTrackHits = 0; 
124   fTrackHitsOld = 0;   
125   fHitType = 2; //default CONTAINERS - based on ROOT structure 
126   fTPCParam = 0;    
127   fNoiseTable = 0;
128   fActiveSectors =0;
129
130 }
131  
132 //_____________________________________________________________________________
133 AliTPC::AliTPC(const char *name, const char *title)
134       : AliDetector(name,title)
135 {
136   //
137   // Standard constructor
138   //
139
140   //
141   // Initialise arrays of hits and digits 
142   fHits     = new TClonesArray("AliTPChit",  176);
143   gAlice->AddHitList(fHits); 
144   fDigitsArray = 0;
145   fClustersArray= 0;
146   fDefaults = 0;
147   //
148   fTrackHits = new AliTPCTrackHitsV2;  
149   fTrackHits->SetHitPrecision(0.002);
150   fTrackHits->SetStepPrecision(0.003);  
151   fTrackHits->SetMaxDistance(100);
152
153   fTrackHitsOld = new AliTPCTrackHits;  //MI - 13.09.2000
154   fTrackHitsOld->SetHitPrecision(0.002);
155   fTrackHitsOld->SetStepPrecision(0.003);  
156   fTrackHitsOld->SetMaxDistance(100); 
157
158   fNoiseTable =0;
159
160   fHitType = 2;
161   fActiveSectors = 0;
162   //
163   // Initialise counters
164   fNsectors = 0;
165
166   //
167   fIshunt     =  0;
168   //
169   // Initialise color attributes
170   SetMarkerColor(kYellow);
171
172   //
173   //  Set TPC parameters
174   //
175
176
177   if (!strcmp(title,"Default")) {       
178     fTPCParam = new AliTPCParamSR;
179   } else {
180     cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
181     fTPCParam=0;
182   }
183
184 }
185
186 //_____________________________________________________________________________
187 AliTPC::~AliTPC()
188 {
189   //
190   // TPC destructor
191   //
192
193   fIshunt   = 0;
194   delete fHits;
195   delete fDigits;
196   delete fTPCParam;
197   delete fTrackHits; //MI 15.09.2000
198   delete fTrackHitsOld; //MI 10.12.2001
199   if (fNoiseTable) delete [] fNoiseTable;
200
201 }
202
203 //_____________________________________________________________________________
204 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
205 {
206   //
207   // Add a hit to the list
208   //
209   //  TClonesArray &lhits = *fHits;
210   //  new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
211   if (fHitType&1){
212     TClonesArray &lhits = *fHits;
213     new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
214   }
215   if (fHitType>1)
216    AddHit2(track,vol,hits);
217 }
218
219 //_____________________________________________________________________________
220 void AliTPC::BuildGeometry()
221 {
222
223   //
224   // Build TPC ROOT TNode geometry for the event display
225   //
226   TNode *nNode, *nTop;
227   TTUBS *tubs;
228   Int_t i;
229   const int kColorTPC=19;
230   char name[5], title[25];
231   const Double_t kDegrad=TMath::Pi()/180;
232   const Double_t kRaddeg=180./TMath::Pi();
233
234
235   Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
236   Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
237
238   Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
239   Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
240
241   Int_t nLo = fTPCParam->GetNInnerSector()/2;
242   Int_t nHi = fTPCParam->GetNOuterSector()/2;  
243
244   const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
245   const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
246   const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
247   const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);  
248
249
250   const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
251   const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
252
253   Double_t rl,ru;
254   
255
256   //
257   // Get ALICE top node
258   //
259
260   nTop=gAlice->GetGeometry()->GetNode("alice");
261
262   //  inner sectors
263
264   rl = fTPCParam->GetInnerRadiusLow();
265   ru = fTPCParam->GetInnerRadiusUp();
266  
267
268   for(i=0;i<nLo;i++) {
269     sprintf(name,"LS%2.2d",i);
270     name[4]='\0';
271     sprintf(title,"TPC low sector %3d",i);
272     title[24]='\0';
273     
274     tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
275                      kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
276     tubs->SetNumberOfDivisions(1);
277     nTop->cd();
278     nNode = new TNode(name,title,name,0,0,0,"");
279     nNode->SetLineColor(kColorTPC);
280     fNodes->Add(nNode);
281   }
282
283   // Outer sectors
284
285   rl = fTPCParam->GetOuterRadiusLow();
286   ru = fTPCParam->GetOuterRadiusUp();
287
288   for(i=0;i<nHi;i++) {
289     sprintf(name,"US%2.2d",i);
290     name[4]='\0';
291     sprintf(title,"TPC upper sector %d",i);
292     title[24]='\0';
293     tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
294                      khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
295     tubs->SetNumberOfDivisions(1);
296     nTop->cd();
297     nNode = new TNode(name,title,name,0,0,0,"");
298     nNode->SetLineColor(kColorTPC);
299     fNodes->Add(nNode);
300   }
301
302 }    
303
304 //_____________________________________________________________________________
305 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
306 {
307   //
308   // Calculate distance from TPC to mouse on the display
309   // Dummy procedure
310   //
311   return 9999;
312 }
313
314 void AliTPC::Clusters2Tracks() 
315  {
316   //-----------------------------------------------------------------
317   // This is a track finder.
318   //-----------------------------------------------------------------
319   AliTPCtracker tracker(fTPCParam,0,fLoader->GetEventFolder()->GetName());
320   tracker.Clusters2Tracks();
321  }
322
323 //_____________________________________________________________________________
324 void AliTPC::CreateMaterials()
325 {
326   //-----------------------------------------------
327   // Create Materials for for TPC simulations
328   //-----------------------------------------------
329
330   //-----------------------------------------------------------------
331   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
332   //-----------------------------------------------------------------
333
334   Int_t iSXFLD=gAlice->Field()->Integ();
335   Float_t sXMGMX=gAlice->Field()->Max();
336
337   Float_t amat[5]; // atomic numbers
338   Float_t zmat[5]; // z
339   Float_t wmat[5]; // proportions
340
341   Float_t density;
342   Float_t apure[2];
343
344
345   //***************** Gases *************************
346   
347   //-------------------------------------------------
348   // pure gases
349   //-------------------------------------------------
350
351   // Neon
352
353
354   amat[0]= 20.18;
355   zmat[0]= 10.;  
356   density = 0.0009;
357  
358   apure[0]=amat[0];
359
360   AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
361
362   // Argon
363
364   amat[0]= 39.948;
365   zmat[0]= 18.;  
366   density = 0.001782;  
367
368   apure[1]=amat[0];
369
370   AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
371  
372
373   //--------------------------------------------------------------
374   // gases - compounds
375   //--------------------------------------------------------------
376
377   Float_t amol[3];
378
379   // CO2
380
381   amat[0]=12.011;
382   amat[1]=15.9994;
383
384   zmat[0]=6.;
385   zmat[1]=8.;
386
387   wmat[0]=1.;
388   wmat[1]=2.;
389
390   density=0.001977;
391
392   amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
393
394   AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
395   
396   // CF4
397
398   amat[0]=12.011;
399   amat[1]=18.998;
400
401   zmat[0]=6.;
402   zmat[1]=9.;
403  
404   wmat[0]=1.;
405   wmat[1]=4.;
406  
407   density=0.003034;
408
409   amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
410
411   AliMixture(11,"CF4",amat,zmat,density,-2,wmat); 
412
413
414   // CH4
415
416   amat[0]=12.011;
417   amat[1]=1.;
418
419   zmat[0]=6.;
420   zmat[1]=1.;
421
422   wmat[0]=1.;
423   wmat[1]=4.;
424
425   density=0.000717;
426
427   amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
428
429   AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
430
431   //----------------------------------------------------------------
432   // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
433   //----------------------------------------------------------------
434
435   char namate[21]; 
436   density = 0.;
437   Float_t am=0;
438   Int_t nc;
439   Float_t rho,absl,X0,buf[1];
440   Int_t nbuf;
441   Float_t a,z;
442
443   for(nc = 0;nc<fNoComp;nc++)
444     {
445     
446       // retrive material constants
447       
448       gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
449
450       amat[nc] = a;
451       zmat[nc] = z;
452
453       Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
454  
455       am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]); 
456       density += fMixtProp[nc]*rho;  // density of the mixture
457       
458     }
459
460   // mixture proportions by weight!
461
462   for(nc = 0;nc<fNoComp;nc++)
463     {
464
465       Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
466
467       wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ? 
468                  apure[nnc] : amol[nnc])/am;
469
470     } 
471
472   // Drift gases 1 - nonsensitive, 2 - sensitive
473
474   AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
475   AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
476
477
478   // Air
479
480   amat[0] = 14.61;
481   zmat[0] = 7.3;
482   density = 0.001205;
483
484   AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.); 
485
486
487   //----------------------------------------------------------------------
488   //               solid materials
489   //----------------------------------------------------------------------
490
491
492   // Kevlar C14H22O2N2
493
494   amat[0] = 12.011;
495   amat[1] = 1.;
496   amat[2] = 15.999;
497   amat[3] = 14.006;
498
499   zmat[0] = 6.;
500   zmat[1] = 1.;
501   zmat[2] = 8.;
502   zmat[3] = 7.;
503
504   wmat[0] = 14.;
505   wmat[1] = 22.;
506   wmat[2] = 2.;
507   wmat[3] = 2.;
508
509   density = 1.45;
510
511   AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);  
512
513   // NOMEX
514
515   amat[0] = 12.011;
516   amat[1] = 1.;
517   amat[2] = 15.999;
518   amat[3] = 14.006;
519
520   zmat[0] = 6.;
521   zmat[1] = 1.;
522   zmat[2] = 8.;
523   zmat[3] = 7.;
524
525   wmat[0] = 14.;
526   wmat[1] = 22.;
527   wmat[2] = 2.;
528   wmat[3] = 2.;
529
530   density = 0.03;
531
532   
533   AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
534
535   // Makrolon C16H18O3
536
537   amat[0] = 12.011;
538   amat[1] = 1.;
539   amat[2] = 15.999;
540
541   zmat[0] = 6.;
542   zmat[1] = 1.;
543   zmat[2] = 8.;
544
545   wmat[0] = 16.;
546   wmat[1] = 18.;
547   wmat[2] = 3.;
548   
549   density = 1.2;
550
551   AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
552   
553   // Mylar C5H4O2
554
555   amat[0]=12.011;
556   amat[1]=1.;
557   amat[2]=15.9994;
558
559   zmat[0]=6.;
560   zmat[1]=1.;
561   zmat[2]=8.;
562
563   wmat[0]=5.;
564   wmat[1]=4.;
565   wmat[2]=2.; 
566
567   density = 1.39;
568   
569   AliMixture(37, "Mylar",amat,zmat,density,-3,wmat); 
570
571   // SiO2 - used later for the glass fiber
572
573   amat[0]=28.086;
574   amat[1]=15.9994;
575
576   zmat[0]=14.;
577   zmat[1]=8.;
578
579   wmat[0]=1.;
580   wmat[1]=2.;
581
582
583   AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
584
585   // Al
586
587   amat[0] = 26.98;
588   zmat[0] = 13.;
589
590   density = 2.7;
591
592   AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
593
594   // Si
595
596   amat[0] = 28.086;
597   zmat[0] = 14.;
598
599   density = 2.33;
600
601   AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
602
603   // Cu
604
605   amat[0] = 63.546;
606   zmat[0] = 29.;
607
608   density = 8.96;
609
610   AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
611
612   // Tedlar C2H3F
613
614   amat[0] = 12.011;
615   amat[1] = 1.;
616   amat[2] = 18.998;
617
618   zmat[0] = 6.;
619   zmat[1] = 1.;
620   zmat[2] = 9.;
621
622   wmat[0] = 2.;
623   wmat[1] = 3.; 
624   wmat[2] = 1.;
625
626   density = 1.71;
627
628   AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);  
629
630
631   // Plexiglas  C5H8O2
632
633   amat[0]=12.011;
634   amat[1]=1.;
635   amat[2]=15.9994;
636
637   zmat[0]=6.;
638   zmat[1]=1.;
639   zmat[2]=8.;
640
641   wmat[0]=5.;
642   wmat[1]=8.;
643   wmat[2]=2.;
644
645   density=1.18;
646
647   AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
648
649   // Epoxy - C14 H20 O3
650
651   
652   amat[0]=12.011;
653   amat[1]=1.;
654   amat[2]=15.9994;
655
656   zmat[0]=6.;
657   zmat[1]=1.;
658   zmat[2]=8.;
659
660   wmat[0]=14.;
661   wmat[1]=20.;
662   wmat[2]=3.;
663
664   density=1.25;
665
666   AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
667
668   // Carbon
669
670   amat[0]=12.011;
671   zmat[0]=6.;
672   density= 2.265;
673
674   AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
675
676   // get epoxy
677
678   gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,X0,absl,buf,nbuf);
679
680   // Carbon fiber
681
682   wmat[0]=0.644; // by weight!
683   wmat[1]=0.356;
684
685   density=0.5*(1.25+2.265);
686
687   AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
688
689   // get SiO2
690
691   gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf); 
692
693   wmat[0]=0.725; // by weight!
694   wmat[1]=0.275;
695
696   density=1.7;
697
698   AliMixture(39,"G10",amat,zmat,density,2,wmat);
699
700  
701
702
703   //----------------------------------------------------------
704   // tracking media for gases
705   //----------------------------------------------------------
706
707   AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
708   AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
709   AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
710   AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001); 
711
712   //-----------------------------------------------------------  
713   // tracking media for solids
714   //-----------------------------------------------------------
715   
716   AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
717   AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
718   AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
719   AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
720   AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
721   AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
722   AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
723   AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
724   AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
725   AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
726   AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
727   AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
728     
729 }
730
731 void AliTPC::GenerNoise(Int_t tablesize)
732 {
733   //
734   //Generate table with noise
735   //
736   if (fTPCParam==0) {
737     // error message
738     fNoiseDepth=0;
739     return;
740   }
741   if (fNoiseTable)  delete[] fNoiseTable;
742   fNoiseTable = new Float_t[tablesize];
743   fNoiseDepth = tablesize; 
744   fCurrentNoise =0; //!index of the noise in  the noise table 
745   
746   Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
747   for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);      
748 }
749
750 Float_t AliTPC::GetNoise()
751 {
752   // get noise from table
753   //  if ((fCurrentNoise%10)==0) 
754   //  fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
755   if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
756   return fNoiseTable[fCurrentNoise++];
757   //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac()); 
758 }
759
760
761 Bool_t  AliTPC::IsSectorActive(Int_t sec)
762 {
763   //
764   // check if the sector is active
765   if (!fActiveSectors) return kTRUE;
766   else return fActiveSectors[sec]; 
767 }
768
769 void    AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
770 {
771   // activate interesting sectors
772   SetTreeAddress();//just for security
773   if (fActiveSectors) delete [] fActiveSectors;
774   fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
775   for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
776   for (Int_t i=0;i<n;i++) 
777     if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector())  fActiveSectors[sectors[i]]=kTRUE;
778     
779 }
780
781 void    AliTPC::SetActiveSectors(Int_t flag)
782 {
783   //
784   // activate sectors which were hitted by tracks 
785   //loop over tracks
786   SetTreeAddress();//just for security
787   if (fHitType==0) return;  // if Clones hit - not short volume ID information
788   if (fActiveSectors) delete [] fActiveSectors;
789   fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
790   if (flag) {
791     for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
792     return;
793   }
794   for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
795   TBranch * branch=0;
796   if (TreeH() == 0x0)
797    {
798      Fatal("SetActiveSectors","Can not find TreeH in folder");
799      return;
800    }
801   if (fHitType>1) branch = TreeH()->GetBranch("TPC2");
802   else branch = TreeH()->GetBranch("TPC");
803   Stat_t ntracks = TreeH()->GetEntries();
804   // loop over all hits
805   cout<<"\nAliTPC::SetActiveSectors():  Got "<<ntracks<<" tracks\n";
806   
807   for(Int_t track=0;track<ntracks;track++)
808    {
809     ResetHits();
810     //
811     if (fTrackHits && fHitType&4) {
812       TBranch * br1 = TreeH()->GetBranch("fVolumes");
813       TBranch * br2 = TreeH()->GetBranch("fNVolumes");
814       br1->GetEvent(track);
815       br2->GetEvent(track);
816       Int_t *volumes = fTrackHits->GetVolumes();
817       for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
818         fActiveSectors[volumes[j]]=kTRUE;
819     }
820     
821     //
822     if (fTrackHitsOld && fHitType&2) {
823       TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
824       br->GetEvent(track);
825       AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
826       for (UInt_t j=0;j<ar->GetSize();j++){
827         fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
828       } 
829     }    
830   }
831   
832 }  
833
834
835
836
837 void AliTPC::Digits2Clusters(Int_t eventnumber)
838 {
839   //-----------------------------------------------------------------
840   // This is a simple cluster finder.
841   //-----------------------------------------------------------------
842   AliTPCclusterer::Digits2Clusters(fTPCParam, fLoader,eventnumber);
843 }
844
845 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
846 extern Double_t SigmaZ2(Double_t, Double_t);
847 //_____________________________________________________________________________
848 void AliTPC::Hits2Clusters(TFile *of, Int_t eventn)
849 {
850   //--------------------------------------------------------
851   // TPC simple cluster generator from hits
852   // obtained from the TPC Fast Simulator
853   // The point errors are taken from the parametrization
854   //--------------------------------------------------------
855
856   //-----------------------------------------------------------------
857   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
858   //-----------------------------------------------------------------
859   // Adopted to Marian's cluster data structure by I.Belikov, CERN,
860   // Jouri.Belikov@cern.ch
861   //----------------------------------------------------------------
862   
863   /////////////////////////////////////////////////////////////////////////////
864   //
865   //---------------------------------------------------------------------
866   //   ALICE TPC Cluster Parameters
867   //--------------------------------------------------------------------
868        
869   
870
871   // Cluster width in rphi
872   const Float_t kACrphi=0.18322;
873   const Float_t kBCrphi=0.59551e-3;
874   const Float_t kCCrphi=0.60952e-1;
875   // Cluster width in z
876   const Float_t kACz=0.19081;
877   const Float_t kBCz=0.55938e-3;
878   const Float_t kCCz=0.30428;
879
880
881   if (!fLoader) {
882      cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
883      return;
884   }
885
886   //if(fDefaults == 0) SetDefaults();
887
888   Float_t sigmaRphi,sigmaZ,clRphi,clZ;
889   //
890   TParticle *particle; // pointer to a given particle
891   AliTPChit *tpcHit; // pointer to a sigle TPC hit
892   Int_t sector;
893   Int_t ipart;
894   Float_t xyz[5];
895   Float_t pl,pt,tanth,rpad,ratio;
896   Float_t cph,sph;
897   
898   //---------------------------------------------------------------
899   //  Get the access to the tracks 
900   //---------------------------------------------------------------
901   
902   TTree *tH = TreeH();
903   if (tH == 0x0)
904    {
905      Fatal("Hits2Clusters","Can not find TreeH in folder");
906      return;
907    }
908   SetTreeAddress();
909   
910   Stat_t ntracks = tH->GetEntries();
911
912   //Switch to the output file
913   
914   if (fLoader->TreeR() == 0x0) fLoader->MakeTree("R");
915   
916   cout<<"fTPCParam->GetTitle() = "<<fTPCParam->GetTitle()<<endl;
917   
918   AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
919   rl->CdGAFile();
920   fTPCParam->Write(fTPCParam->GetTitle());
921
922   AliTPCClustersArray carray;
923   carray.Setup(fTPCParam);
924   carray.SetClusterType("AliTPCcluster");
925   carray.MakeTree(fLoader->TreeR());
926
927   Int_t nclusters=0; //cluster counter
928   
929   //------------------------------------------------------------
930   // Loop over all sectors (72 sectors for 20 deg
931   // segmentation for both lower and upper sectors)
932   // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
933   // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
934   //
935   // First cluster for sector 0 starts at "0"
936   //------------------------------------------------------------
937    
938   for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
939     //MI change
940     fTPCParam->AdjustCosSin(isec,cph,sph);
941     
942     //------------------------------------------------------------
943     // Loop over tracks
944     //------------------------------------------------------------
945     
946     for(Int_t track=0;track<ntracks;track++){
947       ResetHits();
948       tH->GetEvent(track);
949       //
950       //  Get number of the TPC hits
951       //     
952        tpcHit = (AliTPChit*)FirstHit(-1);
953
954       // Loop over hits
955       //
956        while(tpcHit){
957  
958          if (tpcHit->fQ == 0.) {
959            tpcHit = (AliTPChit*) NextHit();
960            continue; //information about track (I.Belikov)
961          }
962         sector=tpcHit->fSector; // sector number
963
964        if(sector != isec){
965          tpcHit = (AliTPChit*) NextHit();
966          continue; 
967        }
968         ipart=tpcHit->Track();
969         particle=gAlice->Particle(ipart);
970         pl=particle->Pz();
971         pt=particle->Pt();
972         if(pt < 1.e-9) pt=1.e-9;
973         tanth=pl/pt;
974         tanth = TMath::Abs(tanth);
975         rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
976         ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
977
978         //   space-point resolutions
979         
980         sigmaRphi=SigmaY2(rpad,tanth,pt);
981         sigmaZ   =SigmaZ2(rpad,tanth   );
982         
983         //   cluster widths
984         
985         clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
986         clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
987         
988         // temporary protection
989         
990         if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
991         if(sigmaZ < 0.) sigmaZ=0.4e-3;
992         if(clRphi < 0.) clRphi=2.5e-3;
993         if(clZ < 0.) clZ=2.5e-5;
994         
995         //
996         
997         //
998         // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
999         // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1000         //
1001         Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
1002         Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
1003         xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi));   // y
1004           Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
1005           fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1006           Float_t ymax=xprim*TMath::Tan(0.5*alpha);
1007           if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim; 
1008         xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
1009           if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z(); 
1010         xyz[2]=sigmaRphi;                                     // fSigmaY2
1011         xyz[3]=sigmaZ;                                        // fSigmaZ2
1012         xyz[4]=tpcHit->fQ;                                    // q
1013
1014         AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
1015         if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);     
1016
1017         Int_t tracks[3]={tpcHit->Track(), -1, -1};
1018         AliTPCcluster cluster(tracks,xyz);
1019
1020         clrow->InsertCluster(&cluster); nclusters++;
1021
1022         tpcHit = (AliTPChit*)NextHit();
1023         
1024
1025       } // end of loop over hits
1026
1027     }   // end of loop over tracks
1028
1029     Int_t nrows=fTPCParam->GetNRow(isec);
1030     for (Int_t irow=0; irow<nrows; irow++) {
1031         AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
1032         if (!clrow) continue;
1033         carray.StoreRow(isec,irow);
1034         carray.ClearRow(isec,irow);
1035     }
1036
1037   } // end of loop over sectors  
1038
1039   cerr<<"Number of made clusters : "<<nclusters<<"                        \n";
1040   fLoader->WriteRecPoints("OVERWRITE");
1041   
1042   
1043 } // end of function
1044
1045 //_________________________________________________________________
1046 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1047 {
1048   //--------------------------------------------------------
1049   //calculate exact cross point of track and given pad row
1050   //resulting values are expressed in "digit" coordinata
1051   //--------------------------------------------------------
1052
1053   //-----------------------------------------------------------------
1054   // Origin: Marian Ivanov  GSI Darmstadt, m.ivanov@gsi.de
1055   //-----------------------------------------------------------------
1056   //
1057   if (fClustersArray==0){    
1058     return;
1059   }
1060   //
1061   TParticle *particle; // pointer to a given particle
1062   AliTPChit *tpcHit; // pointer to a sigle TPC hit
1063   //  Int_t sector,nhits;
1064   Int_t ipart;
1065   const Int_t kcmaxhits=30000;
1066   AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
1067   AliTPCFastVector & xxx = *xxxx;
1068   Int_t maxhits = kcmaxhits;
1069   //construct array for each padrow
1070   for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++) 
1071     fClustersArray->CreateRow(isec,i);
1072   
1073   //---------------------------------------------------------------
1074   //  Get the access to the tracks 
1075   //---------------------------------------------------------------
1076   
1077   TTree *tH = TreeH();
1078   if (tH == 0x0)
1079    {
1080      Fatal("Hits2Clusters","Can not find TreeH in folder");
1081      return;
1082    }
1083   SetTreeAddress();
1084
1085   Stat_t ntracks = tH->GetEntries();
1086   Int_t npart = gAlice->GetNtrack();
1087   //MI change
1088   TBranch * branch=0;
1089   if (fHitType>1) branch = tH->GetBranch("TPC2");
1090   else branch = tH->GetBranch("TPC");
1091
1092   //------------------------------------------------------------
1093   // Loop over tracks
1094   //------------------------------------------------------------
1095
1096   for(Int_t track=0;track<ntracks;track++){ 
1097     Bool_t isInSector=kTRUE;
1098     ResetHits();
1099      isInSector = TrackInVolume(isec,track);
1100     if (!isInSector) continue;
1101     //MI change
1102     branch->GetEntry(track); // get next track
1103     //
1104     //  Get number of the TPC hits and a pointer
1105     //  to the particles
1106     // Loop over hits
1107     //
1108     Int_t currentIndex=0;
1109     Int_t lastrow=-1;  //last writen row
1110
1111     //M.I. changes
1112
1113     tpcHit = (AliTPChit*)FirstHit(-1);
1114     while(tpcHit){
1115       
1116       Int_t sector=tpcHit->fSector; // sector number
1117       if(sector != isec){
1118         tpcHit = (AliTPChit*) NextHit();
1119         continue; 
1120       }
1121
1122       ipart=tpcHit->Track();
1123       if (ipart<npart) particle=gAlice->Particle(ipart);
1124       
1125       //find row number
1126
1127       Float_t  x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1128       Int_t    index[3]={1,isec,0};
1129       Int_t    currentrow = fTPCParam->GetPadRow(x,index) ;     
1130       if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
1131       if (lastrow<0) lastrow=currentrow;
1132       if (currentrow==lastrow){
1133         if ( currentIndex>=maxhits){
1134           maxhits+=kcmaxhits;
1135           xxx.ResizeTo(4*maxhits);
1136         }     
1137         xxx(currentIndex*4)=x[0];
1138         xxx(currentIndex*4+1)=x[1];
1139         xxx(currentIndex*4+2)=x[2];     
1140         xxx(currentIndex*4+3)=tpcHit->fQ;
1141         currentIndex++; 
1142       }
1143       else 
1144         if (currentIndex>2){
1145           Float_t sumx=0;
1146           Float_t sumx2=0;
1147           Float_t sumx3=0;
1148           Float_t sumx4=0;
1149           Float_t sumy=0;
1150           Float_t sumxy=0;
1151           Float_t sumx2y=0;
1152           Float_t sumz=0;
1153           Float_t sumxz=0;
1154           Float_t sumx2z=0;
1155           Float_t sumq=0;
1156           for (Int_t index=0;index<currentIndex;index++){
1157             Float_t x,x2,x3,x4;
1158             x=x2=x3=x4=xxx(index*4);
1159             x2*=x;
1160             x3*=x2;
1161             x4*=x3;
1162             sumx+=x;
1163             sumx2+=x2;
1164             sumx3+=x3;
1165             sumx4+=x4;
1166             sumy+=xxx(index*4+1);
1167             sumxy+=xxx(index*4+1)*x;
1168             sumx2y+=xxx(index*4+1)*x2;
1169             sumz+=xxx(index*4+2);
1170             sumxz+=xxx(index*4+2)*x;
1171             sumx2z+=xxx(index*4+2)*x2;   
1172             sumq+=xxx(index*4+3);
1173           }
1174           Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1175           Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1176             sumx2*(sumx*sumx3-sumx2*sumx2);
1177           
1178           Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1179             sumx2*(sumxy*sumx3-sumx2y*sumx2);
1180           Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1181             sumx2*(sumxz*sumx3-sumx2z*sumx2);
1182           
1183           Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1184             sumx2*(sumx*sumx2y-sumx2*sumxy);
1185           Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1186             sumx2*(sumx*sumx2z-sumx2*sumxz);
1187           
1188           if (TMath::Abs(det)<0.00001){
1189              tpcHit = (AliTPChit*)NextHit();
1190             continue;
1191           }
1192         
1193           Float_t y=detay/det+centralPad;
1194           Float_t z=detaz/det;  
1195           Float_t by=detby/det; //y angle
1196           Float_t bz=detbz/det; //z angle
1197           sumy/=Float_t(currentIndex);
1198           sumz/=Float_t(currentIndex);
1199
1200           AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1201           if (row!=0) {
1202             AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
1203             //    AliComplexCluster cl;
1204             cl->fX=z;
1205             cl->fY=y;
1206             cl->fQ=sumq;
1207             cl->fSigmaX2=bz;
1208             cl->fSigmaY2=by;
1209             cl->fTracks[0]=ipart;
1210           }
1211           currentIndex=0;
1212           lastrow=currentrow;
1213         } //end of calculating cluster for given row
1214         
1215         
1216       tpcHit = (AliTPChit*)NextHit();
1217     } // end of loop over hits
1218   }   // end of loop over tracks 
1219   //write padrows to tree 
1220   for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1221     fClustersArray->StoreRow(isec,ii);    
1222     fClustersArray->ClearRow(isec,ii);        
1223   }
1224   xxxx->Delete();
1225  
1226 }
1227
1228
1229
1230 //__
1231 void AliTPC::SDigits2Digits2(Int_t eventnumber)  
1232 {
1233   //create digits from summable digits
1234   GenerNoise(500000); //create teble with noise
1235
1236   //conect tree with sSDigits
1237   TTree *t = fLoader->TreeS();
1238
1239   if (t == 0x0) 
1240    {
1241      fLoader->LoadSDigits("READ");
1242      t = fLoader->TreeS();
1243      if (t == 0x0)
1244       {
1245         Error("SDigits2Digits2","Can not get input TreeS");
1246         return;
1247       }
1248    }
1249   
1250   if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
1251   
1252   AliSimDigits digarr, *dummy=&digarr;
1253   TBranch* sdb = t->GetBranch("Segment");
1254   if (sdb == 0x0)
1255    {
1256      Error("SDigits2Digits2","Can not find branch with segments in TreeS.");
1257      return;
1258    }  
1259
1260   sdb->SetAddress(&dummy);
1261       
1262   Stat_t nentries = t->GetEntries();
1263
1264   // set zero suppression
1265
1266   fTPCParam->SetZeroSup(2);
1267
1268   // get zero suppression
1269
1270   Int_t zerosup = fTPCParam->GetZeroSup();
1271
1272   //make tree with digits 
1273   
1274   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
1275   arr->SetClass("AliSimDigits");
1276   arr->Setup(fTPCParam);
1277   arr->MakeTree(fLoader->TreeD());
1278   
1279   AliTPCParam * par = fTPCParam;
1280
1281   //Loop over segments of the TPC
1282
1283   for (Int_t n=0; n<nentries; n++) {
1284     t->GetEvent(n);
1285     Int_t sec, row;
1286     if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1287       cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
1288       continue;
1289     }
1290     if (!IsSectorActive(sec)) 
1291      {
1292 //       cout<<n<<" NOT Active \n";
1293        continue;
1294      }
1295     else
1296      {
1297 //       cout<<n<<" Active \n";
1298      }
1299     AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1300     Int_t nrows = digrow->GetNRows();
1301     Int_t ncols = digrow->GetNCols();
1302
1303     digrow->ExpandBuffer();
1304     digarr.ExpandBuffer();
1305     digrow->ExpandTrackBuffer();
1306     digarr.ExpandTrackBuffer();
1307
1308     
1309     Short_t * pamp0 = digarr.GetDigits();
1310     Int_t   * ptracks0 = digarr.GetTracks();
1311     Short_t * pamp1 = digrow->GetDigits();
1312     Int_t   * ptracks1 = digrow->GetTracks();
1313     Int_t  nelems =nrows*ncols;
1314     Int_t saturation = fTPCParam->GetADCSat();
1315     //use internal structure of the AliDigits - for speed reason
1316     //if you cahnge implementation
1317     //of the Alidigits - it must be rewriten -
1318     for (Int_t i= 0; i<nelems; i++){
1319       //      Float_t q = *pamp0;
1320       //q/=16.;  //conversion faktor
1321       //Float_t noise= GetNoise(); 
1322       //q+=noise;      
1323       //q= TMath::Nint(q);
1324       Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1325       if (q>zerosup){
1326         if (q>saturation) q=saturation;      
1327         *pamp1=(Short_t)q;
1328         //if (ptracks0[0]==0)
1329         //  ptracks1[0]=1;
1330         //else
1331         ptracks1[0]=ptracks0[0];        
1332         ptracks1[nelems]=ptracks0[nelems];
1333         ptracks1[2*nelems]=ptracks0[2*nelems];
1334       }
1335       pamp0++;
1336       pamp1++;
1337       ptracks0++;
1338       ptracks1++;        
1339     }
1340
1341     arr->StoreRow(sec,row);
1342     arr->ClearRow(sec,row);   
1343     // cerr<<sec<<"\t"<<row<<"\n";   
1344   }  
1345
1346     
1347   //write results
1348   fLoader->WriteDigits("OVERWRITE");
1349    
1350   delete arr;
1351 }
1352 //__________________________________________________________________
1353 void AliTPC::SetDefaults(){
1354
1355    
1356    cerr<<"Setting default parameters...\n";
1357
1358   // Set response functions
1359
1360   //
1361   AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1362   rl->CdGAFile();
1363   AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1364   if(param){
1365     printf("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...\n");
1366     delete param;
1367     param = new AliTPCParamSR();
1368   }
1369   else {
1370     param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1371   }
1372   if(!param){
1373     printf("No TPC parameters found\n");
1374     exit(4);
1375   }
1376
1377
1378   AliTPCPRF2D    * prfinner   = new AliTPCPRF2D;
1379   AliTPCPRF2D    * prfouter1   = new AliTPCPRF2D;
1380   AliTPCPRF2D    * prfouter2   = new AliTPCPRF2D;  
1381   AliTPCRF1D     * rf    = new AliTPCRF1D(kTRUE);
1382   rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1383   rf->SetOffset(3*param->GetZSigma());
1384   rf->Update();
1385   
1386   TDirectory *savedir=gDirectory;
1387   TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1388   if (!f->IsOpen()) { 
1389     cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1390      exit(3);
1391   }
1392   prfinner->Read("prf_07504_Gati_056068_d02");
1393   prfouter1->Read("prf_10006_Gati_047051_d03");
1394   prfouter2->Read("prf_15006_Gati_047051_d03");  
1395   f->Close();
1396   savedir->cd();
1397
1398   param->SetInnerPRF(prfinner);
1399   param->SetOuter1PRF(prfouter1); 
1400   param->SetOuter2PRF(prfouter2);
1401   param->SetTimeRF(rf);
1402
1403   // set fTPCParam
1404
1405   SetParam(param);
1406
1407
1408   fDefaults = 1;
1409
1410 }
1411 //__________________________________________________________________  
1412 void AliTPC::Hits2Digits(Int_t eventnumber)  
1413
1414  //----------------------------------------------------
1415  // Loop over all sectors for a single event
1416  //----------------------------------------------------
1417   AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1418   rl->GetEvent(eventnumber);
1419   if (fLoader->TreeH() == 0x0)
1420    {
1421      if(fLoader->LoadHits())
1422       {
1423         Error("Hits2Digits","Can not load hits.");
1424       }
1425    }
1426   SetTreeAddress();
1427   
1428   if (fLoader->TreeD() == 0x0 ) 
1429    {
1430      fLoader->MakeTree("D");
1431      if (fLoader->TreeD() == 0x0 ) 
1432       {
1433        Error("Hits2Digits","Can not get TreeD");
1434        return;
1435       }
1436    }
1437
1438   if(fDefaults == 0) SetDefaults();  // check if the parameters are set
1439   GenerNoise(500000); //create teble with noise
1440
1441   //setup TPCDigitsArray 
1442
1443   if(GetDigitsArray()) delete GetDigitsArray();
1444
1445   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
1446   arr->SetClass("AliSimDigits");
1447   arr->Setup(fTPCParam);
1448
1449   arr->MakeTree(fLoader->TreeD());
1450   SetDigitsArray(arr);
1451
1452   fDigitsSwitch=0; // standard digits
1453
1454   cerr<<"Digitizing TPC -- normal digits...\n";
1455
1456  for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) 
1457   if (IsSectorActive(isec)) 
1458    {
1459     cout<<"Sector "<<isec<<"is active\n";
1460     Hits2DigitsSector(isec);
1461    }
1462   else
1463    {
1464     cout<<"Sector "<<isec<<"is NOT active\n";
1465    }
1466
1467   fLoader->WriteDigits("OVERWRITE"); 
1468 }
1469
1470
1471
1472 //__________________________________________________________________
1473 void AliTPC::Hits2SDigits2(Int_t eventnumber)  
1474
1475
1476   //-----------------------------------------------------------
1477   //   summable digits - 16 bit "ADC", no noise, no saturation
1478   //-----------------------------------------------------------
1479
1480  //----------------------------------------------------
1481  // Loop over all sectors for a single event
1482  //----------------------------------------------------
1483 //  AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1484
1485   AliRunLoader* rl = fLoader->GetRunLoader();
1486
1487   rl->GetEvent(eventnumber);
1488   if (fLoader->TreeH() == 0x0)
1489    {
1490      if(fLoader->LoadHits())
1491       {
1492         Error("Hits2Digits","Can not load hits.");
1493         return;
1494       }
1495    }
1496   SetTreeAddress();
1497
1498
1499   if (fLoader->TreeS() == 0x0 ) 
1500    {
1501      fLoader->MakeTree("S");
1502    }
1503   
1504   if(fDefaults == 0) SetDefaults();
1505   
1506   GenerNoise(500000); //create table with noise
1507   //setup TPCDigitsArray 
1508
1509   if(GetDigitsArray()) delete GetDigitsArray();
1510
1511   
1512   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
1513   arr->SetClass("AliSimDigits");
1514   arr->Setup(fTPCParam);
1515   arr->MakeTree(fLoader->TreeS());
1516
1517   SetDigitsArray(arr);
1518
1519   cerr<<"Digitizing TPC -- summable digits...\n"; 
1520
1521   fDigitsSwitch=1; // summable digits
1522   
1523     // set zero suppression to "0"
1524
1525   fTPCParam->SetZeroSup(0);
1526
1527  for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) 
1528   if (IsSectorActive(isec)) 
1529    {
1530 //    cout<<"Sector "<<isec<<" is active\n";
1531     Hits2DigitsSector(isec);
1532    }
1533
1534  fLoader->WriteSDigits("OVERWRITE");
1535
1536 //this line prevents the crash in the similar one
1537 //on the beginning of this method
1538 //destructor attempts to reset the tree, which is deleted by the loader
1539 //need to be redesign
1540  if(GetDigitsArray()) delete GetDigitsArray();
1541  SetDigitsArray(0x0);
1542 }
1543 //__________________________________________________________________
1544
1545 void AliTPC::Hits2SDigits()  
1546
1547
1548   //-----------------------------------------------------------
1549   //   summable digits - 16 bit "ADC", no noise, no saturation
1550   //-----------------------------------------------------------
1551
1552  //----------------------------------------------------
1553  // Loop over all sectors for a single event
1554  //----------------------------------------------------
1555   //MI change - for pp run
1556 //  Int_t eventnumber = gAlice->GetEvNumber();
1557 //  AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1558 //  rl->GetEvent(eventnumber);
1559 // 12/05/2003 This method is obsolete and not used. It should be redesingned
1560 // M.Kowalski
1561
1562   if (fLoader->TreeH() == 0x0)
1563    {
1564      if(fLoader->LoadHits())
1565       {
1566         Error("Hits2Digits","Can not load hits.");
1567       }
1568    }
1569   SetTreeAddress();
1570
1571   if(fDefaults == 0) SetDefaults();
1572   GenerNoise(500000); //create table with noise
1573
1574   //setup TPCDigitsArray 
1575
1576   if(GetDigitsArray()) delete GetDigitsArray();
1577
1578   if (fLoader->TreeS() == 0x0 ) 
1579    {
1580      fLoader->MakeTree("S");
1581    }
1582   
1583   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
1584   arr->SetClass("AliSimDigits");
1585   arr->Setup(fTPCParam);
1586   arr->MakeTree(fLoader->TreeS());
1587   SetDigitsArray(arr);
1588
1589   cerr<<"Digitizing TPC -- summable digits...\n"; 
1590
1591   //  fDigitsSwitch=1; // summable digits  -for the moment direct
1592
1593   for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1594
1595   // write results
1596   //
1597   cout<<"Why method TPC::Hits2SDigits writes digits and not sdigits? skowron\n";
1598   fLoader->WriteDigits("OVERWRITE");
1599 }
1600 //_____________________________________________________________________________
1601
1602 void AliTPC::Hits2DigitsSector(Int_t isec)
1603 {
1604   //-------------------------------------------------------------------
1605   // TPC conversion from hits to digits.
1606   //------------------------------------------------------------------- 
1607
1608   //-----------------------------------------------------------------
1609   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1610   //-----------------------------------------------------------------
1611
1612   //-------------------------------------------------------
1613   //  Get the access to the track hits
1614   //-------------------------------------------------------
1615
1616   // check if the parameters are set - important if one calls this method
1617   // directly, not from the Hits2Digits
1618
1619   if(fDefaults == 0) SetDefaults();
1620
1621   TTree *tH = TreeH(); // pointer to the hits tree
1622   if (tH == 0x0)
1623    {
1624      Fatal("Hits2DigitsSector","Can not find TreeH in folder");
1625      return;
1626    }
1627
1628   Stat_t ntracks = tH->GetEntries();
1629
1630   if( ntracks > 0){
1631
1632   //------------------------------------------- 
1633   //  Only if there are any tracks...
1634   //-------------------------------------------
1635
1636     TObjArray **row;
1637     
1638     //printf("*** Processing sector number %d ***\n",isec);
1639
1640       Int_t nrows =fTPCParam->GetNRow(isec);
1641
1642       row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1643     
1644       MakeSector(isec,nrows,tH,ntracks,row);
1645
1646       //--------------------------------------------------------
1647       //   Digitize this sector, row by row
1648       //   row[i] is the pointer to the TObjArray of AliTPCFastVectors,
1649       //   each one containing electrons accepted on this
1650       //   row, assigned into tracks
1651       //--------------------------------------------------------
1652
1653       Int_t i;
1654
1655       if (fDigitsArray->GetTree()==0) 
1656        {
1657          Fatal("Hits2DigitsSector","Tree not set in fDigitsArray");
1658        }
1659
1660       for (i=0;i<nrows;i++){
1661
1662         AliDigits * dig = fDigitsArray->CreateRow(isec,i); 
1663
1664         DigitizeRow(i,isec,row);
1665
1666         fDigitsArray->StoreRow(isec,i);
1667
1668         Int_t ndig = dig->GetDigitSize(); 
1669         
1670         if (gDebug > 10) 
1671         printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);        
1672         
1673         fDigitsArray->ClearRow(isec,i);  
1674
1675    
1676        } // end of the sector digitization
1677
1678       for(i=0;i<nrows+2;i++){
1679         row[i]->Delete();  
1680         delete row[i];   
1681       }
1682       
1683        delete [] row; // delete the array of pointers to TObjArray-s
1684         
1685   } // ntracks >0
1686
1687 } // end of Hits2DigitsSector
1688
1689
1690 //_____________________________________________________________________________
1691 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1692 {
1693   //-----------------------------------------------------------
1694   // Single row digitization, coupling from the neighbouring
1695   // rows taken into account
1696   //-----------------------------------------------------------
1697
1698   //-----------------------------------------------------------------
1699   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1700   // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1701   //-----------------------------------------------------------------
1702  
1703
1704   Float_t zerosup = fTPCParam->GetZeroSup();
1705   //  Int_t nrows =fTPCParam->GetNRow(isec);
1706   fCurrentIndex[1]= isec;
1707   
1708
1709   Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1710   Int_t nofTbins = fTPCParam->GetMaxTBin();
1711   Int_t indexRange[4];
1712   //
1713   //  Integrated signal for this row
1714   //  and a single track signal
1715   //    
1716
1717   AliTPCFastMatrix *m1 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // integrated
1718   AliTPCFastMatrix *m2 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // single
1719   //
1720   AliTPCFastMatrix &total  = *m1;
1721
1722   //  Array of pointers to the label-signal list
1723
1724   Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1725   Float_t  **pList = new Float_t* [nofDigits]; 
1726
1727   Int_t lp;
1728   Int_t i1;   
1729   for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1730   //
1731   //calculate signal 
1732   //
1733   //Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1734   //Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1735   Int_t row1=irow;
1736   Int_t row2=irow+2; 
1737   for (Int_t row= row1;row<=row2;row++){
1738     Int_t nTracks= rows[row]->GetEntries();
1739     for (i1=0;i1<nTracks;i1++){
1740       fCurrentIndex[2]= row;
1741       fCurrentIndex[3]=irow+1;
1742       if (row==irow+1){
1743         m2->Zero();  // clear single track signal matrix
1744         Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange); 
1745         GetList(trackLabel,nofPads,m2,indexRange,pList);
1746       }
1747       else   GetSignal(rows[row],i1,0,m1,indexRange);
1748     }
1749   }
1750          
1751   Int_t tracks[3];
1752
1753   AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1754   Int_t gi=-1;
1755   Float_t fzerosup = zerosup+0.5;
1756   for(Int_t it=0;it<nofTbins;it++){
1757     Float_t *pq = &(total.UncheckedAt(0,it));
1758     for(Int_t ip=0;ip<nofPads;ip++){
1759       gi++;
1760       Float_t q=*pq;      
1761       pq++;
1762       if(fDigitsSwitch == 0){
1763         q+=GetNoise();
1764         if(q <=fzerosup) continue; // do not fill zeros
1765         q = TMath::Nint(q);
1766         if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat();  // saturation
1767
1768       }
1769
1770       else {
1771        if(q <= 0.) continue; // do not fill zeros
1772        if(q>2000.) q=2000.;
1773        q *= 16.;
1774        q = TMath::Nint(q);
1775       }
1776
1777       //
1778       //  "real" signal or electronic noise (list = -1)?
1779       //    
1780
1781       for(Int_t j1=0;j1<3;j1++){
1782         tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1783       }
1784
1785 //Begin_Html
1786 /*
1787   <A NAME="AliDigits"></A>
1788   using of AliDigits object
1789 */
1790 //End_Html
1791       dig->SetDigitFast((Short_t)q,it,ip);
1792       if (fDigitsArray->IsSimulated())
1793         {
1794          ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1795          ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1796          ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1797         }
1798      
1799     
1800     } // end of loop over time buckets
1801   }  // end of lop over pads 
1802
1803   //
1804   //  This row has been digitized, delete nonused stuff
1805   //
1806
1807   for(lp=0;lp<nofDigits;lp++){
1808     if(pList[lp]) delete [] pList[lp];
1809   }
1810   
1811   delete [] pList;
1812
1813   delete m1;
1814   delete m2;
1815   //  delete m3;
1816
1817 } // end of DigitizeRow
1818
1819 //_____________________________________________________________________________
1820
1821 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, 
1822              AliTPCFastMatrix *m1, AliTPCFastMatrix *m2,Int_t *indexRange)
1823 {
1824
1825   //---------------------------------------------------------------
1826   //  Calculates 2-D signal (pad,time) for a single track,
1827   //  returns a pointer to the signal matrix and the track label 
1828   //  No digitization is performed at this level!!!
1829   //---------------------------------------------------------------
1830
1831   //-----------------------------------------------------------------
1832   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1833   // Modified: Marian Ivanov 
1834   //-----------------------------------------------------------------
1835
1836   AliTPCFastVector *tv;
1837
1838   tv = (AliTPCFastVector*)p1->At(ntr); // pointer to a track
1839   AliTPCFastVector &v = *tv;
1840   
1841   Float_t label = v(0);
1842   Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1)-1)/2;
1843
1844   Int_t nElectrons = (tv->GetNrows()-1)/4;
1845   indexRange[0]=9999; // min pad
1846   indexRange[1]=-1; // max pad
1847   indexRange[2]=9999; //min time
1848   indexRange[3]=-1; // max time
1849
1850   AliTPCFastMatrix &signal = *m1;
1851   AliTPCFastMatrix &total = *m2;
1852   //
1853   //  Loop over all electrons
1854   //
1855   for(Int_t nel=0; nel<nElectrons; nel++){
1856     Int_t idx=nel*4;
1857     Float_t aval =  v(idx+4);
1858     Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac(); 
1859     Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1860     Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
1861
1862     Int_t *index = fTPCParam->GetResBin(0);  
1863     Float_t *weight = & (fTPCParam->GetResWeight(0));
1864
1865     if (n>0) for (Int_t i =0; i<n; i++){       
1866        Int_t pad=index[1]+centralPad;  //in digit coordinates central pad has coordinate 0
1867
1868          if (pad>=0){
1869          Int_t time=index[2];    
1870          Float_t qweight = *(weight)*eltoadcfac;
1871          
1872          if (m1!=0) signal.UncheckedAt(pad,time)+=qweight;
1873          total.UncheckedAt(pad,time)+=qweight;
1874          if (indexRange[0]>pad) indexRange[0]=pad;
1875          if (indexRange[1]<pad) indexRange[1]=pad;
1876          if (indexRange[2]>time) indexRange[2]=time;
1877          if (indexRange[3]<time) indexRange[3]=time;
1878
1879          index+=3;
1880          weight++;      
1881
1882        }         
1883     }
1884   } // end of loop over electrons
1885   
1886   return label; // returns track label when finished
1887 }
1888
1889 //_____________________________________________________________________________
1890 void AliTPC::GetList(Float_t label,Int_t np,AliTPCFastMatrix *m,
1891                      Int_t *indexRange, Float_t **pList)
1892 {
1893   //----------------------------------------------------------------------
1894   //  Updates the list of tracks contributing to digits for a given row
1895   //----------------------------------------------------------------------
1896
1897   //-----------------------------------------------------------------
1898   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1899   //-----------------------------------------------------------------
1900
1901   AliTPCFastMatrix &signal = *m;
1902
1903   // lop over nonzero digits
1904
1905   for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1906     for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1907
1908
1909         // accept only the contribution larger than 500 electrons (1/2 s_noise)
1910
1911         if(signal(ip,it)<0.5) continue; 
1912
1913
1914         Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1915         
1916         if(!pList[globalIndex]){
1917         
1918           // 
1919           // Create new list (6 elements - 3 signals and 3 labels),
1920           //
1921
1922           pList[globalIndex] = new Float_t [6];
1923
1924           // set list to -1 
1925
1926           *pList[globalIndex] = -1.;
1927           *(pList[globalIndex]+1) = -1.;
1928           *(pList[globalIndex]+2) = -1.;
1929           *(pList[globalIndex]+3) = -1.;
1930           *(pList[globalIndex]+4) = -1.;
1931           *(pList[globalIndex]+5) = -1.;
1932
1933
1934           *pList[globalIndex] = label;
1935           *(pList[globalIndex]+3) = signal(ip,it);
1936         }
1937         else{
1938
1939           // check the signal magnitude
1940
1941           Float_t highest = *(pList[globalIndex]+3);
1942           Float_t middle = *(pList[globalIndex]+4);
1943           Float_t lowest = *(pList[globalIndex]+5);
1944
1945           //
1946           //  compare the new signal with already existing list
1947           //
1948
1949           if(signal(ip,it)<lowest) continue; // neglect this track
1950
1951           //
1952
1953           if (signal(ip,it)>highest){
1954             *(pList[globalIndex]+5) = middle;
1955             *(pList[globalIndex]+4) = highest;
1956             *(pList[globalIndex]+3) = signal(ip,it);
1957
1958             *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1959             *(pList[globalIndex]+1) = *pList[globalIndex];
1960             *pList[globalIndex] = label;
1961           }
1962           else if (signal(ip,it)>middle){
1963             *(pList[globalIndex]+5) = middle;
1964             *(pList[globalIndex]+4) = signal(ip,it);
1965
1966             *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1967             *(pList[globalIndex]+1) = label;
1968           }
1969           else{
1970             *(pList[globalIndex]+5) = signal(ip,it);
1971             *(pList[globalIndex]+2) = label;
1972           }
1973         }
1974
1975     } // end of loop over pads
1976   } // end of loop over time bins
1977
1978
1979
1980 }//end of GetList
1981 //___________________________________________________________________
1982 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1983                         Stat_t ntracks,TObjArray **row)
1984 {
1985
1986   //-----------------------------------------------------------------
1987   // Prepares the sector digitization, creates the vectors of
1988   // tracks for each row of this sector. The track vector
1989   // contains the track label and the position of electrons.
1990   //-----------------------------------------------------------------
1991
1992   //-----------------------------------------------------------------
1993   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1994   //-----------------------------------------------------------------
1995
1996   Float_t gasgain = fTPCParam->GetGasGain();
1997   Int_t i;
1998   Float_t xyz[4]; 
1999
2000   AliTPChit *tpcHit; // pointer to a sigle TPC hit    
2001   //MI change
2002   TBranch * branch=0;
2003   if (fHitType>1) branch = TH->GetBranch("TPC2");
2004   else branch = TH->GetBranch("TPC");
2005
2006  
2007   //----------------------------------------------
2008   // Create TObjArray-s, one for each row,
2009   // each TObjArray will store the AliTPCFastVectors
2010   // of electrons, one AliTPCFastVectors per each track.
2011   //---------------------------------------------- 
2012     
2013   Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
2014   AliTPCFastVector **tracks = new AliTPCFastVector* [nrows+2]; //pointers to the track vectors
2015
2016   for(i=0; i<nrows+2; i++){
2017     row[i] = new TObjArray;
2018     nofElectrons[i]=0;
2019     tracks[i]=0;
2020   }
2021
2022  
2023
2024   //--------------------------------------------------------------------
2025   //  Loop over tracks, the "track" contains the full history
2026   //--------------------------------------------------------------------
2027
2028   Int_t previousTrack,currentTrack;
2029   previousTrack = -1; // nothing to store so far!
2030
2031   for(Int_t track=0;track<ntracks;track++){
2032     Bool_t isInSector=kTRUE;
2033     ResetHits();
2034     isInSector = TrackInVolume(isec,track);
2035     if (!isInSector) continue;
2036     //MI change
2037     branch->GetEntry(track); // get next track
2038
2039     //M.I. changes
2040
2041     tpcHit = (AliTPChit*)FirstHit(-1);
2042
2043     //--------------------------------------------------------------
2044     //  Loop over hits
2045     //--------------------------------------------------------------
2046
2047
2048     while(tpcHit){
2049       
2050       Int_t sector=tpcHit->fSector; // sector number
2051       if(sector != isec){
2052         tpcHit = (AliTPChit*) NextHit();
2053         continue; 
2054       }
2055
2056         currentTrack = tpcHit->Track(); // track number
2057
2058
2059         if(currentTrack != previousTrack){
2060                           
2061            // store already filled fTrack
2062               
2063            for(i=0;i<nrows+2;i++){
2064              if(previousTrack != -1){
2065                if(nofElectrons[i]>0){
2066                  AliTPCFastVector &v = *tracks[i];
2067                  v(0) = previousTrack;
2068                  tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2069                  row[i]->Add(tracks[i]);                     
2070                }
2071                else{
2072                  delete tracks[i]; // delete empty AliTPCFastVector
2073                  tracks[i]=0;
2074                }
2075              }
2076
2077              nofElectrons[i]=0;
2078              tracks[i] = new AliTPCFastVector(481); // AliTPCFastVectors for the next fTrack
2079
2080            } // end of loop over rows
2081                
2082            previousTrack=currentTrack; // update track label 
2083         }
2084            
2085         Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2086
2087        //---------------------------------------------------
2088        //  Calculate the electron attachment probability
2089        //---------------------------------------------------
2090
2091
2092         Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
2093                                                         /fTPCParam->GetDriftV(); 
2094         // in microseconds!     
2095         Float_t attProb = fTPCParam->GetAttCoef()*
2096           fTPCParam->GetOxyCont()*time; //  fraction! 
2097    
2098         //-----------------------------------------------
2099         //  Loop over electrons
2100         //-----------------------------------------------
2101         Int_t index[3];
2102         index[1]=isec;
2103         for(Int_t nel=0;nel<qI;nel++){
2104           // skip if electron lost due to the attachment
2105           if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2106           xyz[0]=tpcHit->X();
2107           xyz[1]=tpcHit->Y();
2108           xyz[2]=tpcHit->Z();   
2109           //
2110           // protection for the nonphysical avalanche size (10**6 maximum)
2111           //  
2112           Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2113           xyz[3]= (Float_t) (-gasgain*TMath::Log(rn)); 
2114           index[0]=1;
2115           
2116           TransportElectron(xyz,index);    
2117           Int_t rowNumber;
2118           fTPCParam->GetPadRow(xyz,index); 
2119           // row 0 - cross talk from the innermost row
2120           // row fNRow+1 cross talk from the outermost row
2121           rowNumber = index[2]+1; 
2122           //transform position to local digit coordinates
2123           //relative to nearest pad row 
2124           if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2125           Float_t x1,y1;
2126           if (isec <fTPCParam->GetNInnerSector()) {
2127             x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2128             y1 = fTPCParam->GetYInner(rowNumber);
2129           }
2130           else{
2131             x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2132             y1 = fTPCParam->GetYOuter(rowNumber);
2133           }
2134           // gain inefficiency at the wires edges - linear
2135           x1=TMath::Abs(x1);
2136           y1-=1.;
2137           if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));       
2138        
2139           nofElectrons[rowNumber]++;      
2140           //----------------------------------
2141           // Expand vector if necessary
2142           //----------------------------------
2143           if(nofElectrons[rowNumber]>120){
2144             Int_t range = tracks[rowNumber]->GetNrows();
2145             if((nofElectrons[rowNumber])>(range-1)/4){
2146         
2147               tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
2148             }
2149           }
2150           
2151           AliTPCFastVector &v = *tracks[rowNumber];
2152           Int_t idx = 4*nofElectrons[rowNumber]-3;
2153           Real_t * position = &(((AliTPCFastVector&)v).UncheckedAt(idx)); //make code faster
2154           memcpy(position,xyz,4*sizeof(Float_t));
2155  
2156         } // end of loop over electrons
2157
2158         tpcHit = (AliTPChit*)NextHit();
2159         
2160       } // end of loop over hits
2161     } // end of loop over tracks
2162
2163     //
2164     //   store remaining track (the last one) if not empty
2165     //
2166
2167      for(i=0;i<nrows+2;i++){
2168        if(nofElectrons[i]>0){
2169           AliTPCFastVector &v = *tracks[i];
2170           v(0) = previousTrack;
2171           tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2172           row[i]->Add(tracks[i]);  
2173         }
2174         else{
2175           delete tracks[i];
2176           tracks[i]=0;
2177         }  
2178       }  
2179
2180           delete [] tracks;
2181           delete [] nofElectrons;
2182  
2183
2184 } // end of MakeSector
2185
2186
2187 //_____________________________________________________________________________
2188 void AliTPC::Init()
2189 {
2190   //
2191   // Initialise TPC detector after definition of geometry
2192   //
2193   Int_t i;
2194   //
2195   if(fDebug) {
2196     printf("\n%s: ",ClassName());
2197     for(i=0;i<35;i++) printf("*");
2198     printf(" TPC_INIT ");
2199     for(i=0;i<35;i++) printf("*");
2200     printf("\n%s: ",ClassName());
2201     //
2202     for(i=0;i<80;i++) printf("*");
2203     printf("\n");
2204   }
2205 }
2206
2207 //_____________________________________________________________________________
2208 void AliTPC::MakeBranch(Option_t* option)
2209 {
2210   //
2211   // Create Tree branches for the TPC.
2212   //
2213   if(GetDebug()) Info("MakeBranch","");
2214   Int_t buffersize = 4000;
2215   char branchname[10];
2216   sprintf(branchname,"%s",GetName());
2217   
2218   const char *h = strstr(option,"H");
2219
2220   if ( h && (fHitType<=1) && (fHits == 0x0)) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2221   
2222   AliDetector::MakeBranch(option);
2223
2224   const char *d = strstr(option,"D");
2225  
2226   if (fDigits   && fLoader->TreeD() && d) 
2227    {
2228       MakeBranchInTree(gAlice->TreeD(), branchname, &fDigits, buffersize, 0);
2229    }    
2230
2231   if (fHitType>1) MakeBranch2(option,0); // MI change 14.09.2000
2232 }
2233  
2234 //_____________________________________________________________________________
2235 void AliTPC::ResetDigits()
2236 {
2237   //
2238   // Reset number of digits and the digits array for this detector
2239   //
2240   fNdigits   = 0;
2241   if (fDigits)   fDigits->Clear();
2242 }
2243
2244 //_____________________________________________________________________________
2245 void AliTPC::SetSecAL(Int_t sec)
2246 {
2247   //---------------------------------------------------
2248   // Activate/deactivate selection for lower sectors
2249   //---------------------------------------------------
2250
2251   //-----------------------------------------------------------------
2252   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2253   //-----------------------------------------------------------------
2254   fSecAL = sec;
2255 }
2256
2257 //_____________________________________________________________________________
2258 void AliTPC::SetSecAU(Int_t sec)
2259 {
2260   //----------------------------------------------------
2261   // Activate/deactivate selection for upper sectors
2262   //---------------------------------------------------
2263
2264   //-----------------------------------------------------------------
2265   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2266   //-----------------------------------------------------------------
2267   fSecAU = sec;
2268 }
2269
2270 //_____________________________________________________________________________
2271 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2272 {
2273   //----------------------------------------
2274   // Select active lower sectors
2275   //----------------------------------------
2276
2277   //-----------------------------------------------------------------
2278   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2279   //-----------------------------------------------------------------
2280
2281   fSecLows[0] = s1;
2282   fSecLows[1] = s2;
2283   fSecLows[2] = s3;
2284   fSecLows[3] = s4;
2285   fSecLows[4] = s5;
2286   fSecLows[5] = s6;
2287 }
2288
2289 //_____________________________________________________________________________
2290 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2291                        Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10, 
2292                        Int_t s11 , Int_t s12)
2293 {
2294   //--------------------------------
2295   // Select active upper sectors
2296   //--------------------------------
2297
2298   //-----------------------------------------------------------------
2299   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2300   //-----------------------------------------------------------------
2301
2302   fSecUps[0] = s1;
2303   fSecUps[1] = s2;
2304   fSecUps[2] = s3;
2305   fSecUps[3] = s4;
2306   fSecUps[4] = s5;
2307   fSecUps[5] = s6;
2308   fSecUps[6] = s7;
2309   fSecUps[7] = s8;
2310   fSecUps[8] = s9;
2311   fSecUps[9] = s10;
2312   fSecUps[10] = s11;
2313   fSecUps[11] = s12;
2314 }
2315
2316 //_____________________________________________________________________________
2317 void AliTPC::SetSens(Int_t sens)
2318 {
2319
2320   //-------------------------------------------------------------
2321   // Activates/deactivates the sensitive strips at the center of
2322   // the pad row -- this is for the space-point resolution calculations
2323   //-------------------------------------------------------------
2324
2325   //-----------------------------------------------------------------
2326   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2327   //-----------------------------------------------------------------
2328
2329   fSens = sens;
2330 }
2331
2332  
2333 void AliTPC::SetSide(Float_t side=0.)
2334 {
2335   // choice of the TPC side
2336
2337   fSide = side;
2338  
2339 }
2340 //____________________________________________________________________________
2341 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2342                            Float_t p2,Float_t p3)
2343 {
2344
2345   // gax mixture definition
2346
2347  fNoComp = nc;
2348  
2349  fMixtComp[0]=c1;
2350  fMixtComp[1]=c2;
2351  fMixtComp[2]=c3;
2352
2353  fMixtProp[0]=p1;
2354  fMixtProp[1]=p2;
2355  fMixtProp[2]=p3; 
2356  
2357  
2358 }
2359 //_____________________________________________________________________________
2360
2361 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2362 {
2363   //
2364   // electron transport taking into account:
2365   // 1. diffusion, 
2366   // 2.ExB at the wires
2367   // 3. nonisochronity
2368   //
2369   // xyz and index must be already transformed to system 1
2370   //
2371
2372   fTPCParam->Transform1to2(xyz,index);
2373   
2374   //add diffusion
2375   Float_t driftl=xyz[2];
2376   if(driftl<0.01) driftl=0.01;
2377   driftl=TMath::Sqrt(driftl);
2378   Float_t sigT = driftl*(fTPCParam->GetDiffT());
2379   Float_t sigL = driftl*(fTPCParam->GetDiffL());
2380   xyz[0]=gRandom->Gaus(xyz[0],sigT);
2381   xyz[1]=gRandom->Gaus(xyz[1],sigT);
2382   xyz[2]=gRandom->Gaus(xyz[2],sigL);
2383
2384   // ExB
2385   
2386   if (fTPCParam->GetMWPCReadout()==kTRUE){
2387     Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2388     xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2389   }
2390   //add nonisochronity (not implemented yet)  
2391 }
2392   
2393 ClassImp(AliTPCdigit)
2394  
2395 //_____________________________________________________________________________
2396 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2397   AliDigit(tracks)
2398 {
2399   //
2400   // Creates a TPC digit object
2401   //
2402   fSector     = digits[0];
2403   fPadRow     = digits[1];
2404   fPad        = digits[2];
2405   fTime       = digits[3];
2406   fSignal     = digits[4];
2407 }
2408
2409  
2410 ClassImp(AliTPChit)
2411  
2412 //_____________________________________________________________________________
2413 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2414 AliHit(shunt,track)
2415 {
2416   //
2417   // Creates a TPC hit object
2418   //
2419   fSector     = vol[0];
2420   fPadRow     = vol[1];
2421   fX          = hits[0];
2422   fY          = hits[1];
2423   fZ          = hits[2];
2424   fQ          = hits[3];
2425 }
2426  
2427 //________________________________________________________________________
2428 // Additional code because of the AliTPCTrackHitsV2
2429
2430 void AliTPC::MakeBranch2(Option_t *option,const char *file)
2431 {
2432   //
2433   // Create a new branch in the current Root Tree
2434   // The branch of fHits is automatically split
2435   // MI change 14.09.2000
2436   if(GetDebug()) Info("MakeBranch2","");
2437   if (fHitType<2) return;
2438   char branchname[10];
2439   sprintf(branchname,"%s2",GetName());  
2440   //
2441   // Get the pointer to the header
2442   const char *cH = strstr(option,"H");
2443   //
2444   if (fTrackHits   && TreeH() && cH && fHitType&4) 
2445    {
2446     if(GetDebug()) Info("MakeBranch2","Making branch for Type 4 Hits");
2447     TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2448    }    
2449
2450   if (fTrackHitsOld   && TreeH() && cH && fHitType&2) 
2451    {    
2452     if(GetDebug()) Info("MakeBranch2","Making branch for Type 2 Hits");
2453     AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld, 
2454                                                    TreeH(),fBufferSize,99);
2455     TreeH()->GetListOfBranches()->Add(branch);
2456    }    
2457 }
2458
2459 void AliTPC::SetTreeAddress()
2460 {
2461 //Sets tree address for hits  
2462   if (fHitType<=1)
2463    {
2464      if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2465      AliDetector::SetTreeAddress();
2466    }
2467   if (fHitType>1) SetTreeAddress2();
2468 }
2469
2470 void AliTPC::SetTreeAddress2()
2471 {
2472   //
2473   // Set branch address for the TrackHits Tree
2474   // 
2475   if(GetDebug()) Info("SetTreeAddress2","");
2476   
2477   TBranch *branch;
2478   char branchname[20];
2479   sprintf(branchname,"%s2",GetName());
2480   //
2481   // Branch address for hit tree
2482   TTree *treeH = TreeH();
2483   if ((treeH)&&(fHitType&4)) {
2484     branch = treeH->GetBranch(branchname);
2485     if (branch) 
2486      {
2487        branch->SetAddress(&fTrackHits);
2488        cout<<"AliTPC::SetTreeAddress2 fHitType&4 Setting"<<endl;       
2489      }
2490     else cout<<"AliTPC::SetTreeAddress2 fHitType&4  Failed"<<endl;
2491   }
2492   if ((treeH)&&(fHitType&2)) {
2493     branch = treeH->GetBranch(branchname);
2494     if (branch) 
2495      {
2496        branch->SetAddress(&fTrackHitsOld);
2497        cout<<"AliTPC::SetTreeAddress2 fHitType&2 Setting"<<endl;       
2498      }
2499     else cout<<"AliTPC::SetTreeAddress2 fHitType&2  Failed"<<endl;
2500   }
2501   //set address to TREETR
2502
2503   TTree *treeTR = TreeTR();
2504   if (treeTR && fTrackReferences) {
2505     branch = treeTR->GetBranch(GetName());
2506     if (branch) branch->SetAddress(&fTrackReferences);
2507   }
2508
2509 }
2510
2511 void AliTPC::FinishPrimary()
2512 {
2513   if (fTrackHits &&fHitType&4)      fTrackHits->FlushHitStack();  
2514   if (fTrackHitsOld && fHitType&2)  fTrackHitsOld->FlushHitStack();  
2515 }
2516
2517
2518 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2519
2520   //
2521   // add hit to the list  
2522   Int_t rtrack;
2523   if (fIshunt) {
2524     int primary = gAlice->GetPrimary(track);
2525     gAlice->Particle(primary)->SetBit(kKeepBit);
2526     rtrack=primary;
2527   } else {
2528     rtrack=track;
2529     gAlice->FlagTrack(track);
2530   }  
2531   //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2532   //if (hit->fTrack!=rtrack)
2533   //  cout<<"bad track number\n";
2534   if (fTrackHits && fHitType&4) 
2535     fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2536                              hits[1],hits[2],(Int_t)hits[3]);
2537   if (fTrackHitsOld &&fHitType&2 ) 
2538     fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2539                                 hits[1],hits[2],(Int_t)hits[3]);
2540   
2541 }
2542
2543 void AliTPC::ResetHits()
2544
2545   if (fHitType&1) AliDetector::ResetHits();
2546   if (fHitType>1) ResetHits2();
2547 }
2548
2549 void AliTPC::ResetHits2()
2550 {
2551   //
2552   //reset hits
2553   if (fTrackHits && fHitType&4) fTrackHits->Clear();
2554   if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2555
2556 }   
2557
2558 AliHit* AliTPC::FirstHit(Int_t track)
2559 {
2560   if (fHitType>1) return FirstHit2(track);
2561   return AliDetector::FirstHit(track);
2562 }
2563 AliHit* AliTPC::NextHit()
2564 {
2565   if (fHitType>1) return NextHit2();
2566   
2567   return AliDetector::NextHit();
2568 }
2569
2570 AliHit* AliTPC::FirstHit2(Int_t track)
2571 {
2572   //
2573   // Initialise the hit iterator
2574   // Return the address of the first hit for track
2575   // If track>=0 the track is read from disk
2576   // while if track<0 the first hit of the current
2577   // track is returned
2578   // 
2579   if(track>=0) {
2580     gAlice->ResetHits();
2581     TreeH()->GetEvent(track);
2582   }
2583   //
2584   if (fTrackHits && fHitType&4) {
2585     fTrackHits->First();
2586     return fTrackHits->GetHit();
2587   }
2588   if (fTrackHitsOld && fHitType&2) {
2589     fTrackHitsOld->First();
2590     return fTrackHitsOld->GetHit();
2591   }
2592
2593   else return 0;
2594 }
2595
2596 AliHit* AliTPC::NextHit2()
2597 {
2598   //
2599   //Return the next hit for the current track
2600
2601
2602   if (fTrackHitsOld && fHitType&2) {
2603     fTrackHitsOld->Next();
2604     return fTrackHitsOld->GetHit();
2605   }
2606   if (fTrackHits) {
2607     fTrackHits->Next();
2608     return fTrackHits->GetHit();
2609   }
2610   else 
2611     return 0;
2612 }
2613
2614 void AliTPC::LoadPoints(Int_t)
2615 {
2616   //
2617   Int_t a = 0;
2618   /*  if(fHitType==1) return AliDetector::LoadPoints(a);
2619   LoadPoints2(a);
2620   */
2621   if(fHitType==1) AliDetector::LoadPoints(a);
2622   else LoadPoints2(a);
2623    
2624   // LoadPoints3(a);
2625
2626 }
2627
2628
2629 void AliTPC::RemapTrackHitIDs(Int_t *map)
2630 {
2631   if (!fTrackHits) return;
2632   
2633   if (fTrackHitsOld && fHitType&2){
2634     AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2635     for (UInt_t i=0;i<arr->GetSize();i++){
2636       AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2637       info->fTrackID = map[info->fTrackID];
2638     }
2639   }
2640   if (fTrackHitsOld && fHitType&4){
2641     TClonesArray * arr = fTrackHits->GetArray();;
2642     for (Int_t i=0;i<arr->GetEntriesFast();i++){
2643       AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2644       info->fTrackID = map[info->fTrackID];
2645     }
2646   }
2647 }
2648
2649 Bool_t   AliTPC::TrackInVolume(Int_t id,Int_t track)
2650 {
2651   //return bool information - is track in given volume
2652   //load only part of the track information 
2653   //return true if current track is in volume
2654   //
2655   //  return kTRUE;
2656   if (fTrackHitsOld && fHitType&2) {
2657     TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
2658     br->GetEvent(track);
2659     AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2660     for (UInt_t j=0;j<ar->GetSize();j++){
2661       if (  ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2662     } 
2663   }
2664
2665   if (fTrackHits && fHitType&4) {
2666     TBranch * br1 = TreeH()->GetBranch("fVolumes");
2667     TBranch * br2 = TreeH()->GetBranch("fNVolumes");    
2668     br2->GetEvent(track);
2669     br1->GetEvent(track);    
2670     Int_t *volumes = fTrackHits->GetVolumes();
2671     Int_t nvolumes = fTrackHits->GetNVolumes();
2672     if (!volumes && nvolumes>0) {
2673       printf("Problematic track\t%d\t%d",track,nvolumes);
2674       return kFALSE;
2675     }
2676     for (Int_t j=0;j<nvolumes; j++)
2677       if (volumes[j]==id) return kTRUE;    
2678   }
2679
2680   if (fHitType&1) {
2681     TBranch * br = TreeH()->GetBranch("fSector");
2682     br->GetEvent(track);
2683     for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2684       if (  ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2685     } 
2686   }
2687   return kFALSE;  
2688
2689 }
2690
2691 //_____________________________________________________________________________
2692 void AliTPC::LoadPoints2(Int_t)
2693 {
2694   //
2695   // Store x, y, z of all hits in memory
2696   //
2697   if (fTrackHits == 0 && fTrackHitsOld==0) return;
2698   //
2699   Int_t nhits =0;
2700   if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2701   if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2702   
2703   if (nhits == 0) return;
2704   Int_t tracks = gAlice->GetNtrack();
2705   if (fPoints == 0) fPoints = new TObjArray(tracks);
2706   AliHit *ahit;
2707   //
2708   Int_t *ntrk=new Int_t[tracks];
2709   Int_t *limi=new Int_t[tracks];
2710   Float_t **coor=new Float_t*[tracks];
2711   for(Int_t i=0;i<tracks;i++) {
2712     ntrk[i]=0;
2713     coor[i]=0;
2714     limi[i]=0;
2715   }
2716   //
2717   AliPoints *points = 0;
2718   Float_t *fp=0;
2719   Int_t trk;
2720   Int_t chunk=nhits/4+1;
2721   //
2722   // Loop over all the hits and store their position
2723   //
2724   ahit = FirstHit2(-1);
2725   while (ahit){
2726     trk=ahit->GetTrack();
2727     if(ntrk[trk]==limi[trk]) {
2728       //
2729       // Initialise a new track
2730       fp=new Float_t[3*(limi[trk]+chunk)];
2731       if(coor[trk]) {
2732         memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2733         delete [] coor[trk];
2734       }
2735       limi[trk]+=chunk;
2736       coor[trk] = fp;
2737     } else {
2738       fp = coor[trk];
2739     }
2740     fp[3*ntrk[trk]  ] = ahit->X();
2741     fp[3*ntrk[trk]+1] = ahit->Y();
2742     fp[3*ntrk[trk]+2] = ahit->Z();
2743     ntrk[trk]++;
2744     ahit = NextHit2();
2745   }
2746
2747
2748
2749   //
2750   for(trk=0; trk<tracks; ++trk) {
2751     if(ntrk[trk]) {
2752       points = new AliPoints();
2753       points->SetMarkerColor(GetMarkerColor());
2754       points->SetMarkerSize(GetMarkerSize());
2755       points->SetDetector(this);
2756       points->SetParticle(trk);
2757       points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2758       fPoints->AddAt(points,trk);
2759       delete [] coor[trk];
2760       coor[trk]=0;
2761     }
2762   }
2763   delete [] coor;
2764   delete [] ntrk;
2765   delete [] limi;
2766 }
2767
2768
2769 //_____________________________________________________________________________
2770 void AliTPC::LoadPoints3(Int_t)
2771 {
2772   //
2773   // Store x, y, z of all hits in memory
2774   // - only intersection point with pad row
2775   if (fTrackHits == 0) return;
2776   //
2777   Int_t nhits = fTrackHits->GetEntriesFast();
2778   if (nhits == 0) return;
2779   Int_t tracks = gAlice->GetNtrack();
2780   if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2781   fPoints->Expand(2*tracks);
2782   AliHit *ahit;
2783   //
2784   Int_t *ntrk=new Int_t[tracks];
2785   Int_t *limi=new Int_t[tracks];
2786   Float_t **coor=new Float_t*[tracks];
2787   for(Int_t i=0;i<tracks;i++) {
2788     ntrk[i]=0;
2789     coor[i]=0;
2790     limi[i]=0;
2791   }
2792   //
2793   AliPoints *points = 0;
2794   Float_t *fp=0;
2795   Int_t trk;
2796   Int_t chunk=nhits/4+1;
2797   //
2798   // Loop over all the hits and store their position
2799   //
2800   ahit = FirstHit2(-1);
2801   //for (Int_t hit=0;hit<nhits;hit++) {
2802
2803   Int_t lastrow = -1;
2804   while (ahit){
2805     //    ahit = (AliHit*)fHits->UncheckedAt(hit);
2806     trk=ahit->GetTrack(); 
2807     Float_t  x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2808     Int_t    index[3]={1,((AliTPChit*)ahit)->fSector,0};
2809     Int_t    currentrow = fTPCParam->GetPadRow(x,index) ;
2810     if (currentrow!=lastrow){
2811       lastrow = currentrow;
2812       //later calculate intersection point           
2813       if(ntrk[trk]==limi[trk]) {
2814         //
2815         // Initialise a new track
2816         fp=new Float_t[3*(limi[trk]+chunk)];
2817         if(coor[trk]) {
2818           memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2819           delete [] coor[trk];
2820         }
2821         limi[trk]+=chunk;
2822         coor[trk] = fp;
2823       } else {
2824         fp = coor[trk];
2825       }
2826       fp[3*ntrk[trk]  ] = ahit->X();
2827       fp[3*ntrk[trk]+1] = ahit->Y();
2828       fp[3*ntrk[trk]+2] = ahit->Z();
2829       ntrk[trk]++;
2830     }
2831     ahit = NextHit2();
2832   }
2833   
2834   //
2835   for(trk=0; trk<tracks; ++trk) {
2836     if(ntrk[trk]) {
2837       points = new AliPoints();
2838       points->SetMarkerColor(GetMarkerColor()+1);
2839       points->SetMarkerStyle(5);
2840       points->SetMarkerSize(0.2);
2841       points->SetDetector(this);
2842       points->SetParticle(trk);
2843       //      points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2844       points->SetPolyMarker(ntrk[trk],coor[trk],30);
2845       fPoints->AddAt(points,tracks+trk);
2846       delete [] coor[trk];
2847       coor[trk]=0;
2848     }
2849   }
2850   delete [] coor;
2851   delete [] ntrk;
2852   delete [] limi;
2853 }
2854
2855
2856
2857 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2858 {
2859
2860   //
2861   //fill clones array with intersection of current point with the
2862   //middle of the row
2863   Int_t sector;
2864   Int_t ipart;
2865   
2866   const Int_t kcmaxhits=30000;
2867   AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
2868   AliTPCFastVector & xxx = *xxxx;
2869   Int_t maxhits = kcmaxhits;
2870       
2871   //
2872   AliTPChit * tpcHit=0;
2873   tpcHit = (AliTPChit*)FirstHit2(-1);
2874   Int_t currentIndex=0;
2875   Int_t lastrow=-1;  //last writen row
2876
2877   while (tpcHit){
2878     if (tpcHit==0) continue;
2879     sector=tpcHit->fSector; // sector number
2880     ipart=tpcHit->Track();
2881     
2882     //find row number
2883     
2884     Float_t  x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2885     Int_t    index[3]={1,sector,0};
2886     Int_t    currentrow = fTPCParam->GetPadRow(x,index) ;       
2887     if (currentrow<0) continue;
2888     if (lastrow<0) lastrow=currentrow;
2889     if (currentrow==lastrow){
2890       if ( currentIndex>=maxhits){
2891         maxhits+=kcmaxhits;
2892         xxx.ResizeTo(4*maxhits);
2893       }     
2894       xxx(currentIndex*4)=x[0];
2895       xxx(currentIndex*4+1)=x[1];
2896       xxx(currentIndex*4+2)=x[2];       
2897       xxx(currentIndex*4+3)=tpcHit->fQ;
2898       currentIndex++;   
2899     }
2900     else 
2901       if (currentIndex>2){
2902         Float_t sumx=0;
2903         Float_t sumx2=0;
2904         Float_t sumx3=0;
2905         Float_t sumx4=0;
2906         Float_t sumy=0;
2907         Float_t sumxy=0;
2908         Float_t sumx2y=0;
2909         Float_t sumz=0;
2910         Float_t sumxz=0;
2911         Float_t sumx2z=0;
2912         Float_t sumq=0;
2913         for (Int_t index=0;index<currentIndex;index++){
2914           Float_t x,x2,x3,x4;
2915           x=x2=x3=x4=xxx(index*4);
2916           x2*=x;
2917           x3*=x2;
2918           x4*=x3;
2919           sumx+=x;
2920           sumx2+=x2;
2921           sumx3+=x3;
2922           sumx4+=x4;
2923           sumy+=xxx(index*4+1);
2924           sumxy+=xxx(index*4+1)*x;
2925           sumx2y+=xxx(index*4+1)*x2;
2926           sumz+=xxx(index*4+2);
2927           sumxz+=xxx(index*4+2)*x;
2928           sumx2z+=xxx(index*4+2)*x2;     
2929           sumq+=xxx(index*4+3);
2930         }
2931         Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2932         Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2933           sumx2*(sumx*sumx3-sumx2*sumx2);
2934         
2935         Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2936           sumx2*(sumxy*sumx3-sumx2y*sumx2);
2937         Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2938           sumx2*(sumxz*sumx3-sumx2z*sumx2);
2939         
2940         Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2941           sumx2*(sumx*sumx2y-sumx2*sumxy);
2942         Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2943           sumx2*(sumx*sumx2z-sumx2*sumxz);
2944         
2945         Float_t y=detay/det+centralPad;
2946         Float_t z=detaz/det;    
2947         Float_t by=detby/det; //y angle
2948         Float_t bz=detbz/det; //z angle
2949         sumy/=Float_t(currentIndex);
2950         sumz/=Float_t(currentIndex);
2951         
2952         AliComplexCluster cl;
2953         cl.fX=z;
2954         cl.fY=y;
2955         cl.fQ=sumq;
2956         cl.fSigmaX2=bz;
2957         cl.fSigmaY2=by;
2958         cl.fTracks[0]=ipart;
2959         
2960         AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2961         if (row!=0) row->InsertCluster(&cl);
2962         currentIndex=0;
2963         lastrow=currentrow;
2964       } //end of calculating cluster for given row
2965                 
2966   } // end of loop over hits
2967   xxxx->Delete();
2968
2969 }
2970 //_______________________________________________________________________________
2971 void AliTPC::Digits2Reco(Int_t firstevent,Int_t lastevent)
2972 {
2973   // produces rec points from digits and writes them on the root file
2974   // AliTPCclusters.root
2975
2976   TDirectory *cwd = gDirectory;
2977
2978
2979   AliTPCParamSR *dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60");
2980   if(dig){
2981     printf("You are running 2 pad-length geom hits with 3 pad-length geom digits\n");
2982     delete dig;
2983     dig = new AliTPCParamSR();
2984   }
2985   else
2986   {
2987    dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60"); 
2988   }
2989   if(!dig){
2990    printf("No TPC parameters found\n");
2991    exit(3);
2992   }
2993    
2994   SetParam(dig);
2995   cout<<"AliTPC::Digits2Reco: TPC parameteres have been set"<<endl; 
2996   TFile *out;
2997   if(!gSystem->Getenv("CONFIG_FILE")){
2998     out=TFile::Open("AliTPCclusters.root","recreate");
2999   }
3000   else {
3001     const char *tmp1;
3002     const char *tmp2;
3003     char tmp3[80];
3004     tmp1=gSystem->Getenv("CONFIG_FILE_PREFIX");
3005     tmp2=gSystem->Getenv("CONFIG_OUTDIR");
3006     sprintf(tmp3,"%s%s/AliTPCclusters.root",tmp1,tmp2);
3007     out=TFile::Open(tmp3,"recreate");
3008   }
3009
3010   TStopwatch timer;
3011   cout<<"AliTPC::Digits2Reco - determination of rec points begins"<<endl;
3012   timer.Start();
3013   cwd->cd();
3014   for(Int_t iev=firstevent;iev<lastevent+1;iev++){
3015
3016     printf("Processing event %d\n",iev);
3017     Digits2Clusters(iev);
3018   }
3019   cout<<"AliTPC::Digits2Reco - determination of rec points ended"<<endl;
3020   timer.Stop();
3021   timer.Print();
3022   out->Close();
3023   cwd->cd(); 
3024
3025
3026 }
3027
3028 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
3029 {
3030  cout<<"AliTPC::MakeLoader ";
3031  
3032  fLoader = new AliTPCLoader(GetName(),topfoldername);
3033  
3034  if (fLoader)
3035   {
3036    cout<<"Success"<<endl;
3037   }
3038  else
3039   {
3040    cout<<"Failure"<<endl;
3041   }
3042
3043  return fLoader;
3044 }
3045
3046 ////////////////////////////////////////////////////////////////////////
3047 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
3048 //
3049 // load TPC paarmeters from a given file or create new if the object
3050 // is not found there
3051 // 12/05/2003 This method should be moved to the AliTPCLoader
3052 // and one has to decide where to store the TPC parameters
3053 // M.Kowalski
3054   char paramName[50];
3055   sprintf(paramName,"75x40_100x60_150x60");
3056   AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
3057   if (paramTPC) {
3058     cout<<"TPC parameters "<<paramName<<" found."<<endl;
3059   } else {
3060     cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3061         <<endl;    
3062     paramTPC = new AliTPCParamSR;
3063   }
3064   return paramTPC;
3065
3066 // the older version of parameters can be accessed with this code.
3067 // In some cases, we have old parameters saved in the file but 
3068 // digits were created with new parameters, it can be distinguish 
3069 // by the name of TPC TreeD. The code here is just for the case 
3070 // we would need to compare with old data, uncomment it if needed.
3071 //
3072 //  char paramName[50];
3073 //  sprintf(paramName,"75x40_100x60");
3074 //  AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
3075 //  if (paramTPC) {
3076 //    cout<<"TPC parameters "<<paramName<<" found."<<endl;
3077 //  } else {
3078 //    sprintf(paramName,"75x40_100x60_150x60");
3079 //    paramTPC=(AliTPCParam*)in->Get(paramName);
3080 //    if (paramTPC) {
3081 //      cout<<"TPC parameters "<<paramName<<" found."<<endl;
3082 //    } else {
3083 //      cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3084 //          <<endl;    
3085 //      paramTPC = new AliTPCParamSR;
3086 //    }
3087 //  }
3088 //  return paramTPC;
3089
3090 }
3091