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