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