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