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