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