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