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