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