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