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