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