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