1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.67 2003/03/03 16:36:16 kowal2
19 Added LoadTPCParam method
21 Revision 1.66 2003/02/11 16:54:07 hristov
22 Updated AliTrackReference class (S.Radomski)
24 Revision 1.65 2002/11/21 22:43:32 alibrary
25 Removing AliMC and AliMCProcess
27 Revision 1.64 2002/10/23 07:17:33 alibrary
28 Introducing Riostream.h
30 Revision 1.63 2002/10/14 14:57:42 hristov
31 Merging the VirtualMC branch to the main development branch (HEAD)
33 Revision 1.54.4.3 2002/10/11 08:34:48 hristov
34 Updating VirtualMC to v3-09-02
36 Revision 1.62 2002/09/23 09:22:56 hristov
37 The address of the TrackReferences is set (M.Ivanov)
39 Revision 1.61 2002/09/10 07:06:42 kowal2
40 Corrected for the memory leak. Thanks to J. Chudoba and M. Ivanov
42 Revision 1.60 2002/06/12 14:56:56 kowal2
43 Added track length to the reference hits
45 Revision 1.59 2002/06/05 15:37:31 kowal2
46 Added cross-talk from the wires beyond the first and the last rows
48 Revision 1.58 2002/05/27 14:33:14 hristov
49 The new class AliTrackReference used (M.Ivanov)
51 Revision 1.57 2002/05/07 17:23:11 kowal2
52 Linear gain inefficiency instead of the step one at the wire edges.
54 Revision 1.56 2002/04/04 16:26:33 kowal2
55 Digits (Sdigits) go to separate files now.
57 Revision 1.55 2002/03/29 06:57:45 kowal2
58 Restored backward compatibility to use the hits from Dec. 2000 production.
60 Revision 1.54 2002/03/18 17:59:13 kowal2
61 Chnges in the pad geometry - 3 pad lengths introduced.
63 Revision 1.53 2002/02/25 11:02:56 kowal2
64 Changes towards speeding up the code. Thanks to Marian Ivanov.
66 Revision 1.52 2002/02/18 09:26:09 kowal2
67 Removed compiler warning
69 Revision 1.51 2002/01/21 17:13:21 kowal2
70 New track hits using root containers. Setting active sectors added.
72 Revision 1.50 2001/12/06 14:16:19 kowal2
75 Revision 1.49 2001/11/30 11:55:37 hristov
76 Noise table created in Hits2SDigits (M.Ivanov)
78 Revision 1.48 2001/11/24 16:10:21 kowal2
81 Revision 1.47 2001/11/19 10:25:34 kowal2
82 Nearest integer instead of integer when converting to ADC counts
84 Revision 1.46 2001/11/07 06:47:12 kowal2
87 Revision 1.45 2001/11/03 13:33:48 kowal2
88 Updated algorithms in Hits2SDigits, SDigits2Digits,
92 Revision 1.44 2001/08/30 09:28:48 hristov
93 TTree names are explicitly set via SetName(name) and then Write() is called
95 Revision 1.43 2001/07/28 12:02:54 hristov
96 Branch split level set to 99
98 Revision 1.42 2001/07/28 11:38:52 hristov
99 Loop variable declared once
101 Revision 1.41 2001/07/28 10:53:50 hristov
102 Digitisation done according to the general scheme (M.Ivanov)
104 Revision 1.40 2001/07/27 13:03:14 hristov
105 Default Branch split level set to 99
107 Revision 1.39 2001/07/26 09:09:34 kowal2
108 Hits2Reco method added
110 Revision 1.38 2001/07/20 14:32:43 kowal2
111 Processing of many events possible now
113 Revision 1.37 2001/06/12 07:17:18 kowal2
114 Hits2SDigits method implemented (summable digits)
116 Revision 1.36 2001/05/16 14:57:25 alibrary
117 New files for folders and Stack
119 Revision 1.35 2001/05/08 16:02:22 kowal2
120 Updated material specifications
122 Revision 1.34 2001/05/08 15:00:15 hristov
123 Corrections for tracking in arbitrary magnenetic field. Changes towards a concept of global Alice track. Back propagation of reconstructed tracks (Yu.Belikov)
125 Revision 1.33 2001/04/03 12:40:43 kowal2
128 Revision 1.32 2001/03/12 17:47:36 hristov
129 Changes needed on Sun with CC 5.0
131 Revision 1.31 2001/03/12 08:21:50 kowal2
132 Corrected C++ bug in the material definitions
134 Revision 1.30 2001/03/01 17:34:47 kowal2
135 Correction due to the accuracy problem
137 Revision 1.29 2001/02/28 16:34:40 kowal2
138 Protection against nonphysical values of the avalanche size,
141 Revision 1.28 2001/01/26 19:57:19 hristov
142 Major upgrade of AliRoot code
144 Revision 1.27 2001/01/13 17:29:33 kowal2
145 Sun compiler correction
147 Revision 1.26 2001/01/10 07:59:43 kowal2
148 Corrections to load points from the noncompressed hits.
150 Revision 1.25 2000/11/02 07:25:31 kowal2
151 Changes due to the new hit structure.
154 Revision 1.24 2000/10/05 16:06:09 kowal2
155 Forward declarations. Changes due to a new class AliComplexCluster.
157 Revision 1.23 2000/10/02 21:28:18 fca
158 Removal of useless dependecies via forward declarations
160 Revision 1.22 2000/07/10 20:57:39 hristov
161 Update of TPC code and macros by M.Kowalski
163 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
164 Changes to obey the coding rules
166 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
167 Splitted from AliTPCtracking
169 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
170 Changed parameter settings
172 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
174 Defaults loaded automatically (hard-wired)
175 Optional parameters can be set via macro called in the constructor
177 Revision 1.19 2000/04/18 19:00:59 fca
178 Small bug fixes to TPC files
180 Revision 1.18 2000/04/17 09:37:33 kowal2
181 removed obsolete AliTPCDigitsDisplay.C
183 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
185 New, experimental data structure from M. Ivanov
186 New tracking algorithm
187 Different pad geometry for different sectors
188 Digitization rewritten
190 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
191 Not used anymore - removed
193 Revision 1.17 2000/01/19 17:17:30 fca
194 Introducing a list of lists of hits -- more hits allowed for detector now
196 Revision 1.16 1999/11/05 09:29:23 fca
197 Accept only signals > 0
199 Revision 1.15 1999/10/08 06:26:53 fca
200 Removed ClustersIndex - not used anymore
202 Revision 1.14 1999/09/29 09:24:33 fca
203 Introduction of the Copyright and cvs Log
207 ///////////////////////////////////////////////////////////////////////////////
209 // Time Projection Chamber //
210 // This class contains the basic functions for the Time Projection Chamber //
211 // detector. Functions specific to one particular geometry are //
212 // contained in the derived classes //
216 <img src="picts/AliTPCClass.gif">
221 ///////////////////////////////////////////////////////////////////////////////
229 #include <TGeometry.h>
232 #include <TObjectTable.h>
233 #include "TParticle.h"
239 #include <Riostream.h>
241 #include <Riostream.h>
243 #include "AliTrackReference.h"
246 #include "AliTPCParamSR.h"
247 #include "AliTPCPRF2D.h"
248 #include "AliTPCRF1D.h"
249 #include "AliDigits.h"
250 #include "AliSimDigits.h"
251 #include "AliTPCTrackHits.h"
252 #include "AliTPCTrackHitsV2.h"
253 #include "AliPoints.h"
254 #include "AliArrayBranch.h"
257 #include "AliTPCDigitsArray.h"
258 #include "AliComplexCluster.h"
259 #include "AliClusters.h"
260 #include "AliTPCClustersRow.h"
261 #include "AliTPCClustersArray.h"
263 #include "AliTPCcluster.h"
264 #include "AliTPCclusterer.h"
265 #include "AliTPCtracker.h"
267 #include <TInterpreter.h>
274 //_____________________________________________________________________________
275 // helper class for fast matrix and vector manipulation - no range checking
276 // origin - Marian Ivanov
278 class AliTPCFastMatrix : public TMatrix {
280 AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb);
281 inline Float_t & UncheckedAt(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
282 inline Float_t UncheckedAtFast(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
285 AliTPCFastMatrix::AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb):
286 TMatrix(row_lwb, row_upb,col_lwb,col_upb)
289 //_____________________________________________________________________________
290 class AliTPCFastVector : public TVector {
292 AliTPCFastVector(Int_t size);
293 virtual ~AliTPCFastVector(){;}
294 inline Float_t & UncheckedAt(Int_t index) const {return fElements[index];} //fast acces
298 AliTPCFastVector::AliTPCFastVector(Int_t size):TVector(size){
301 //_____________________________________________________________________________
305 // Default constructor
316 fHitType = 2; //default CONTAINERS - based on ROOT structure
323 //_____________________________________________________________________________
324 AliTPC::AliTPC(const char *name, const char *title)
325 : AliDetector(name,title)
328 // Standard constructor
332 // Initialise arrays of hits and digits
333 fHits = new TClonesArray("AliTPChit", 176);
334 gAlice->AddHitList(fHits);
339 fTrackHits = new AliTPCTrackHitsV2;
340 fTrackHits->SetHitPrecision(0.002);
341 fTrackHits->SetStepPrecision(0.003);
342 fTrackHits->SetMaxDistance(100);
344 fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
345 fTrackHitsOld->SetHitPrecision(0.002);
346 fTrackHitsOld->SetStepPrecision(0.003);
347 fTrackHitsOld->SetMaxDistance(100);
354 // Initialise counters
360 // Initialise color attributes
361 SetMarkerColor(kYellow);
364 // Set TPC parameters
368 if (!strcmp(title,"Default")) {
369 fTPCParam = new AliTPCParamSR;
371 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
377 //_____________________________________________________________________________
388 delete fTrackHits; //MI 15.09.2000
389 delete fTrackHitsOld; //MI 10.12.2001
390 if (fNoiseTable) delete [] fNoiseTable;
394 //_____________________________________________________________________________
395 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
398 // Add a hit to the list
400 // TClonesArray &lhits = *fHits;
401 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
403 TClonesArray &lhits = *fHits;
404 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
407 AddHit2(track,vol,hits);
410 //_____________________________________________________________________________
411 void AliTPC::BuildGeometry()
415 // Build TPC ROOT TNode geometry for the event display
420 const int kColorTPC=19;
421 char name[5], title[25];
422 const Double_t kDegrad=TMath::Pi()/180;
423 const Double_t kRaddeg=180./TMath::Pi();
426 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
427 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
429 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
430 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
432 Int_t nLo = fTPCParam->GetNInnerSector()/2;
433 Int_t nHi = fTPCParam->GetNOuterSector()/2;
435 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
436 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
437 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
438 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
441 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
442 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
448 // Get ALICE top node
451 nTop=gAlice->GetGeometry()->GetNode("alice");
455 rl = fTPCParam->GetInnerRadiusLow();
456 ru = fTPCParam->GetInnerRadiusUp();
460 sprintf(name,"LS%2.2d",i);
462 sprintf(title,"TPC low sector %3d",i);
465 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
466 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
467 tubs->SetNumberOfDivisions(1);
469 nNode = new TNode(name,title,name,0,0,0,"");
470 nNode->SetLineColor(kColorTPC);
476 rl = fTPCParam->GetOuterRadiusLow();
477 ru = fTPCParam->GetOuterRadiusUp();
480 sprintf(name,"US%2.2d",i);
482 sprintf(title,"TPC upper sector %d",i);
484 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
485 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
486 tubs->SetNumberOfDivisions(1);
488 nNode = new TNode(name,title,name,0,0,0,"");
489 nNode->SetLineColor(kColorTPC);
495 //_____________________________________________________________________________
496 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
499 // Calculate distance from TPC to mouse on the display
505 void AliTPC::Clusters2Tracks(TFile *of) {
506 //-----------------------------------------------------------------
507 // This is a track finder.
508 //-----------------------------------------------------------------
509 AliTPCtracker tracker(fTPCParam);
510 tracker.Clusters2Tracks(gFile,of);
513 //_____________________________________________________________________________
514 void AliTPC::CreateMaterials()
516 //-----------------------------------------------
517 // Create Materials for for TPC simulations
518 //-----------------------------------------------
520 //-----------------------------------------------------------------
521 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
522 //-----------------------------------------------------------------
524 Int_t iSXFLD=gAlice->Field()->Integ();
525 Float_t sXMGMX=gAlice->Field()->Max();
527 Float_t amat[5]; // atomic numbers
528 Float_t zmat[5]; // z
529 Float_t wmat[5]; // proportions
535 //***************** Gases *************************
537 //-------------------------------------------------
539 //-------------------------------------------------
550 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
560 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
563 //--------------------------------------------------------------
565 //--------------------------------------------------------------
582 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
584 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
599 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
601 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
617 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
619 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
621 //----------------------------------------------------------------
622 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
623 //----------------------------------------------------------------
629 Float_t rho,absl,X0,buf[1];
633 for(nc = 0;nc<fNoComp;nc++)
636 // retrive material constants
638 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
643 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
645 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
646 density += fMixtProp[nc]*rho; // density of the mixture
650 // mixture proportions by weight!
652 for(nc = 0;nc<fNoComp;nc++)
655 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
657 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
658 apure[nnc] : amol[nnc])/am;
662 // Drift gases 1 - nonsensitive, 2 - sensitive
664 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
665 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
674 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
677 //----------------------------------------------------------------------
679 //----------------------------------------------------------------------
701 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
723 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
741 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
759 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
761 // SiO2 - used later for the glass fiber
773 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
782 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
791 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
800 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
818 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
837 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
839 // Epoxy - C14 H20 O3
856 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
864 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
868 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,X0,absl,buf,nbuf);
872 wmat[0]=0.644; // by weight!
875 density=0.5*(1.25+2.265);
877 AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
881 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
883 wmat[0]=0.725; // by weight!
888 AliMixture(39,"G10",amat,zmat,density,2,wmat);
893 //----------------------------------------------------------
894 // tracking media for gases
895 //----------------------------------------------------------
897 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
898 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
899 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
900 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
902 //-----------------------------------------------------------
903 // tracking media for solids
904 //-----------------------------------------------------------
906 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
907 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
908 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
909 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
910 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
911 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
912 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
913 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
914 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
915 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
916 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
917 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
921 void AliTPC::GenerNoise(Int_t tablesize)
924 //Generate table with noise
931 if (fNoiseTable) delete[] fNoiseTable;
932 fNoiseTable = new Float_t[tablesize];
933 fNoiseDepth = tablesize;
934 fCurrentNoise =0; //!index of the noise in the noise table
936 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
937 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
940 Float_t AliTPC::GetNoise()
942 // get noise from table
943 // if ((fCurrentNoise%10)==0)
944 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
945 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
946 return fNoiseTable[fCurrentNoise++];
947 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
951 Bool_t AliTPC::IsSectorActive(Int_t sec)
954 // check if the sector is active
955 if (!fActiveSectors) return kTRUE;
956 else return fActiveSectors[sec];
959 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
961 // activate interesting sectors
962 if (fActiveSectors) delete [] fActiveSectors;
963 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
964 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
965 for (Int_t i=0;i<n;i++)
966 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
970 void AliTPC::SetActiveSectors(Int_t flag)
973 // activate sectors which were hitted by tracks
975 if (fHitType==0) return; // if Clones hit - not short volume ID information
976 if (fActiveSectors) delete [] fActiveSectors;
977 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
979 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
982 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
984 if (fHitType>1) branch = gAlice->TreeH()->GetBranch("TPC2");
985 else branch = gAlice->TreeH()->GetBranch("TPC");
986 Stat_t ntracks = gAlice->TreeH()->GetEntries();
987 // loop over all hits
988 for(Int_t track=0;track<ntracks;track++){
991 if (fTrackHits && fHitType&4) {
992 TBranch * br1 = gAlice->TreeH()->GetBranch("fVolumes");
993 TBranch * br2 = gAlice->TreeH()->GetBranch("fNVolumes");
994 br1->GetEvent(track);
995 br2->GetEvent(track);
996 Int_t *volumes = fTrackHits->GetVolumes();
997 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
998 fActiveSectors[volumes[j]]=kTRUE;
1002 if (fTrackHitsOld && fHitType&2) {
1003 TBranch * br = gAlice->TreeH()->GetBranch("fTrackHitsInfo");
1004 br->GetEvent(track);
1005 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
1006 for (UInt_t j=0;j<ar->GetSize();j++){
1007 fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
1017 void AliTPC::Digits2Clusters(TFile *of, Int_t eventnumber)
1019 //-----------------------------------------------------------------
1020 // This is a simple cluster finder.
1021 //-----------------------------------------------------------------
1022 AliTPCclusterer::Digits2Clusters(fTPCParam,of,eventnumber);
1025 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
1026 extern Double_t SigmaZ2(Double_t, Double_t);
1027 //_____________________________________________________________________________
1028 void AliTPC::Hits2Clusters(TFile *of, Int_t eventn)
1030 //--------------------------------------------------------
1031 // TPC simple cluster generator from hits
1032 // obtained from the TPC Fast Simulator
1033 // The point errors are taken from the parametrization
1034 //--------------------------------------------------------
1036 //-----------------------------------------------------------------
1037 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1038 //-----------------------------------------------------------------
1039 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
1040 // Jouri.Belikov@cern.ch
1041 //----------------------------------------------------------------
1043 /////////////////////////////////////////////////////////////////////////////
1045 //---------------------------------------------------------------------
1046 // ALICE TPC Cluster Parameters
1047 //--------------------------------------------------------------------
1051 // Cluster width in rphi
1052 const Float_t kACrphi=0.18322;
1053 const Float_t kBCrphi=0.59551e-3;
1054 const Float_t kCCrphi=0.60952e-1;
1055 // Cluster width in z
1056 const Float_t kACz=0.19081;
1057 const Float_t kBCz=0.55938e-3;
1058 const Float_t kCCz=0.30428;
1060 TDirectory *savedir=gDirectory;
1062 if (!of->IsOpen()) {
1063 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
1067 //if(fDefaults == 0) SetDefaults();
1069 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
1071 TParticle *particle; // pointer to a given particle
1072 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1076 Float_t pl,pt,tanth,rpad,ratio;
1079 //---------------------------------------------------------------
1080 // Get the access to the tracks
1081 //---------------------------------------------------------------
1083 TTree *tH = gAlice->TreeH();
1084 Stat_t ntracks = tH->GetEntries();
1086 //Switch to the output file
1091 sprintf(cname,"TreeC_TPC_%d",eventn);
1093 fTPCParam->Write(fTPCParam->GetTitle());
1094 AliTPCClustersArray carray;
1095 carray.Setup(fTPCParam);
1096 carray.SetClusterType("AliTPCcluster");
1099 Int_t nclusters=0; //cluster counter
1101 //------------------------------------------------------------
1102 // Loop over all sectors (72 sectors for 20 deg
1103 // segmentation for both lower and upper sectors)
1104 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
1105 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1107 // First cluster for sector 0 starts at "0"
1108 //------------------------------------------------------------
1110 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
1112 fTPCParam->AdjustCosSin(isec,cph,sph);
1114 //------------------------------------------------------------
1116 //------------------------------------------------------------
1118 for(Int_t track=0;track<ntracks;track++){
1120 tH->GetEvent(track);
1122 // Get number of the TPC hits
1124 tpcHit = (AliTPChit*)FirstHit(-1);
1130 if (tpcHit->fQ == 0.) {
1131 tpcHit = (AliTPChit*) NextHit();
1132 continue; //information about track (I.Belikov)
1134 sector=tpcHit->fSector; // sector number
1137 tpcHit = (AliTPChit*) NextHit();
1140 ipart=tpcHit->Track();
1141 particle=gAlice->Particle(ipart);
1144 if(pt < 1.e-9) pt=1.e-9;
1146 tanth = TMath::Abs(tanth);
1147 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
1148 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
1150 // space-point resolutions
1152 sigmaRphi=SigmaY2(rpad,tanth,pt);
1153 sigmaZ =SigmaZ2(rpad,tanth );
1157 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
1158 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
1160 // temporary protection
1162 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
1163 if(sigmaZ < 0.) sigmaZ=0.4e-3;
1164 if(clRphi < 0.) clRphi=2.5e-3;
1165 if(clZ < 0.) clZ=2.5e-5;
1170 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
1171 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1173 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
1174 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
1175 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
1176 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
1177 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1178 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
1179 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
1180 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
1181 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
1182 xyz[2]=sigmaRphi; // fSigmaY2
1183 xyz[3]=sigmaZ; // fSigmaZ2
1184 xyz[4]=tpcHit->fQ; // q
1186 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
1187 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
1189 Int_t tracks[3]={tpcHit->Track(), -1, -1};
1190 AliTPCcluster cluster(tracks,xyz);
1192 clrow->InsertCluster(&cluster); nclusters++;
1194 tpcHit = (AliTPChit*)NextHit();
1197 } // end of loop over hits
1199 } // end of loop over tracks
1201 Int_t nrows=fTPCParam->GetNRow(isec);
1202 for (Int_t irow=0; irow<nrows; irow++) {
1203 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
1204 if (!clrow) continue;
1205 carray.StoreRow(isec,irow);
1206 carray.ClearRow(isec,irow);
1209 } // end of loop over sectors
1211 cerr<<"Number of made clusters : "<<nclusters<<" \n";
1212 carray.GetTree()->SetName(cname);
1213 carray.GetTree()->Write();
1215 savedir->cd(); //switch back to the input file
1217 } // end of function
1219 //_________________________________________________________________
1220 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1222 //--------------------------------------------------------
1223 //calculate exact cross point of track and given pad row
1224 //resulting values are expressed in "digit" coordinata
1225 //--------------------------------------------------------
1227 //-----------------------------------------------------------------
1228 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1229 //-----------------------------------------------------------------
1231 if (fClustersArray==0){
1235 TParticle *particle; // pointer to a given particle
1236 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1237 // Int_t sector,nhits;
1239 const Int_t kcmaxhits=30000;
1240 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
1241 AliTPCFastVector & xxx = *xxxx;
1242 Int_t maxhits = kcmaxhits;
1243 //construct array for each padrow
1244 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
1245 fClustersArray->CreateRow(isec,i);
1247 //---------------------------------------------------------------
1248 // Get the access to the tracks
1249 //---------------------------------------------------------------
1251 TTree *tH = gAlice->TreeH();
1252 Stat_t ntracks = tH->GetEntries();
1253 Int_t npart = gAlice->GetNtrack();
1256 if (fHitType>1) branch = tH->GetBranch("TPC2");
1257 else branch = tH->GetBranch("TPC");
1259 //------------------------------------------------------------
1261 //------------------------------------------------------------
1263 for(Int_t track=0;track<ntracks;track++){
1264 Bool_t isInSector=kTRUE;
1266 isInSector = TrackInVolume(isec,track);
1267 if (!isInSector) continue;
1269 branch->GetEntry(track); // get next track
1271 // Get number of the TPC hits and a pointer
1275 Int_t currentIndex=0;
1276 Int_t lastrow=-1; //last writen row
1280 tpcHit = (AliTPChit*)FirstHit(-1);
1283 Int_t sector=tpcHit->fSector; // sector number
1285 tpcHit = (AliTPChit*) NextHit();
1289 ipart=tpcHit->Track();
1290 if (ipart<npart) particle=gAlice->Particle(ipart);
1294 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1295 Int_t index[3]={1,isec,0};
1296 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1297 if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
1298 if (lastrow<0) lastrow=currentrow;
1299 if (currentrow==lastrow){
1300 if ( currentIndex>=maxhits){
1302 xxx.ResizeTo(4*maxhits);
1304 xxx(currentIndex*4)=x[0];
1305 xxx(currentIndex*4+1)=x[1];
1306 xxx(currentIndex*4+2)=x[2];
1307 xxx(currentIndex*4+3)=tpcHit->fQ;
1311 if (currentIndex>2){
1323 for (Int_t index=0;index<currentIndex;index++){
1325 x=x2=x3=x4=xxx(index*4);
1333 sumy+=xxx(index*4+1);
1334 sumxy+=xxx(index*4+1)*x;
1335 sumx2y+=xxx(index*4+1)*x2;
1336 sumz+=xxx(index*4+2);
1337 sumxz+=xxx(index*4+2)*x;
1338 sumx2z+=xxx(index*4+2)*x2;
1339 sumq+=xxx(index*4+3);
1341 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1342 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1343 sumx2*(sumx*sumx3-sumx2*sumx2);
1345 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1346 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1347 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1348 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1350 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1351 sumx2*(sumx*sumx2y-sumx2*sumxy);
1352 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1353 sumx2*(sumx*sumx2z-sumx2*sumxz);
1355 if (TMath::Abs(det)<0.00001){
1356 tpcHit = (AliTPChit*)NextHit();
1360 Float_t y=detay/det+centralPad;
1361 Float_t z=detaz/det;
1362 Float_t by=detby/det; //y angle
1363 Float_t bz=detbz/det; //z angle
1364 sumy/=Float_t(currentIndex);
1365 sumz/=Float_t(currentIndex);
1367 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1369 AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
1370 // AliComplexCluster cl;
1376 cl->fTracks[0]=ipart;
1380 } //end of calculating cluster for given row
1383 tpcHit = (AliTPChit*)NextHit();
1384 } // end of loop over hits
1385 } // end of loop over tracks
1386 //write padrows to tree
1387 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1388 fClustersArray->StoreRow(isec,ii);
1389 fClustersArray->ClearRow(isec,ii);
1398 void AliTPC::SDigits2Digits2(Int_t eventnumber)
1400 //create digits from summable digits
1401 GenerNoise(500000); //create teble with noise
1404 sprintf(sname,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1405 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1407 //conect tree with sSDigits
1409 if (gAlice->GetTreeDFile()) {
1410 t = (TTree *)gAlice->GetTreeDFile()->Get(sname);
1412 t = (TTree *)gDirectory->Get(sname);
1415 cerr<<"TPC tree with sdigits not found"<<endl;
1418 AliSimDigits digarr, *dummy=&digarr;
1419 t->GetBranch("Segment")->SetAddress(&dummy);
1420 Stat_t nentries = t->GetEntries();
1422 // set zero suppression
1424 fTPCParam->SetZeroSup(2);
1426 // get zero suppression
1428 Int_t zerosup = fTPCParam->GetZeroSup();
1430 //make tree with digits
1432 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1433 arr->SetClass("AliSimDigits");
1434 arr->Setup(fTPCParam);
1435 // Note that methods arr->MakeTree have different signatures
1436 if (gAlice->GetTreeDFile()) {
1437 arr->MakeTree(gAlice->GetTreeDFile());
1439 arr->MakeTree(fDigitsFile);
1442 AliTPCParam * par =fTPCParam;
1444 //Loop over segments of the TPC
1446 for (Int_t n=0; n<nentries; n++) {
1449 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1450 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
1453 if (!IsSectorActive(sec)) continue;
1454 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1455 Int_t nrows = digrow->GetNRows();
1456 Int_t ncols = digrow->GetNCols();
1458 digrow->ExpandBuffer();
1459 digarr.ExpandBuffer();
1460 digrow->ExpandTrackBuffer();
1461 digarr.ExpandTrackBuffer();
1464 Short_t * pamp0 = digarr.GetDigits();
1465 Int_t * ptracks0 = digarr.GetTracks();
1466 Short_t * pamp1 = digrow->GetDigits();
1467 Int_t * ptracks1 = digrow->GetTracks();
1468 Int_t nelems =nrows*ncols;
1469 Int_t saturation = fTPCParam->GetADCSat();
1470 //use internal structure of the AliDigits - for speed reason
1471 //if you cahnge implementation
1472 //of the Alidigits - it must be rewriten -
1473 for (Int_t i= 0; i<nelems; i++){
1474 // Float_t q = *pamp0;
1475 //q/=16.; //conversion faktor
1476 //Float_t noise= GetNoise();
1478 //q= TMath::Nint(q);
1479 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1481 if (q>saturation) q=saturation;
1483 //if (ptracks0[0]==0)
1486 ptracks1[0]=ptracks0[0];
1487 ptracks1[nelems]=ptracks0[nelems];
1488 ptracks1[2*nelems]=ptracks0[2*nelems];
1496 arr->StoreRow(sec,row);
1497 arr->ClearRow(sec,row);
1498 // cerr<<sec<<"\t"<<row<<"\n";
1505 arr->GetTree()->SetName(dname);
1506 arr->GetTree()->AutoSave();
1509 //_________________________________________
1510 void AliTPC::Merge(TTree * intree, Int_t *mask, Int_t nin, Int_t outid)
1513 //intree - pointer to array of trees with s digits
1514 //mask - mask for each
1515 //nin - number of inputs
1516 //outid - event number of the output event
1518 //create digits from summable digits
1519 //conect tree with sSDigits
1522 AliSimDigits ** digarr = new AliSimDigits*[nin];
1523 for (Int_t i1=0;i1<nin; i1++){
1525 intree[i1].GetBranch("Segment")->SetAddress(&digarr[i1]);
1527 Stat_t nentries = intree[0].GetEntries();
1529 //make tree with digits
1531 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),outid);
1532 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1533 arr->SetClass("AliSimDigits");
1534 arr->Setup(fTPCParam);
1535 arr->MakeTree(fDigitsFile);
1537 // set zero suppression
1539 fTPCParam->SetZeroSup(2);
1541 // get zero suppression
1543 Int_t zerosup = fTPCParam->GetZeroSup();
1546 AliTPCParam * par =fTPCParam;
1548 //Loop over segments of the TPC
1549 for (Int_t n=0; n<nentries; n++) {
1551 for (Int_t i=0;i<nin; i++){
1552 //connect proper digits
1553 intree[i].GetEvent(n);
1554 digarr[i]->ExpandBuffer();
1555 digarr[i]->ExpandTrackBuffer();
1558 if (!par->AdjustSectorRow(digarr[0]->GetID(),sec,row)) {
1559 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr[0]->GetID()<<endl;
1563 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1564 Int_t nrows = digrow->GetNRows();
1565 Int_t ncols = digrow->GetNCols();
1567 digrow->ExpandBuffer();
1568 digrow->ExpandTrackBuffer();
1570 for (Int_t rows=0;rows<nrows; rows++){
1571 for (Int_t col=0;col<ncols; col++){
1573 Int_t label[1000]; // i hope no more than 300 events merged
1575 // looop over digits
1576 for (Int_t i=0;i<nin; i++){
1577 q += digarr[i]->GetDigitFast(rows,col);
1578 for (Int_t tr=0;tr<3;tr++) {
1579 Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
1581 label[labptr]=lab+mask[i]-2;
1587 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac()*16);
1589 q/=16; //conversion faktor
1594 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat();
1595 digrow->SetDigitFast((Short_t)q,rows,col);
1596 for (Int_t tr=0;tr<3;tr++){
1598 ((AliSimDigits*)digrow)->SetTrackIDFast(label[tr],rows,col,tr);
1600 ((AliSimDigits*)digrow)->SetTrackIDFast(-1,rows,col,tr);
1605 arr->StoreRow(sec,row);
1607 arr->ClearRow(sec,row);
1612 arr->GetTree()->SetName(dname);
1613 arr->GetTree()->Write();
1619 /*_________________________________________
1620 void AliTPC::SDigits2Digits(Int_t eventnumber)
1624 cerr<<"Digitizing TPC...\n";
1626 Hits2Digits(eventnumber);
1631 // char treeName[100];
1633 // sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1635 // GetDigitsArray()->GetTree()->Write(treeName);
1638 //__________________________________________________________________
1639 void AliTPC::SetDefaults(){
1642 cerr<<"Setting default parameters...\n";
1644 // Set response functions
1646 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1648 printf("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...\n");
1650 param = new AliTPCParamSR();
1653 param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1656 printf("No TPC parameters found\n");
1661 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1662 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1663 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1664 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1665 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1666 rf->SetOffset(3*param->GetZSigma());
1669 TDirectory *savedir=gDirectory;
1670 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1672 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1675 prfinner->Read("prf_07504_Gati_056068_d02");
1676 prfouter1->Read("prf_10006_Gati_047051_d03");
1677 prfouter2->Read("prf_15006_Gati_047051_d03");
1681 param->SetInnerPRF(prfinner);
1682 param->SetOuter1PRF(prfouter1);
1683 param->SetOuter2PRF(prfouter2);
1684 param->SetTimeRF(rf);
1694 //__________________________________________________________________
1695 void AliTPC::Hits2Digits(Int_t eventnumber)
1697 //----------------------------------------------------
1698 // Loop over all sectors for a single event
1699 //----------------------------------------------------
1702 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1703 GenerNoise(500000); //create teble with noise
1705 //setup TPCDigitsArray
1707 if(GetDigitsArray()) delete GetDigitsArray();
1709 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1710 arr->SetClass("AliSimDigits");
1711 arr->Setup(fTPCParam);
1712 // Note that methods arr->MakeTree have different signatures
1713 if (gAlice->GetTreeDFile()) {
1714 arr->MakeTree(gAlice->GetTreeDFile());
1716 arr->MakeTree(fDigitsFile);
1718 SetDigitsArray(arr);
1720 fDigitsSwitch=0; // standard digits
1722 cerr<<"Digitizing TPC -- normal digits...\n";
1724 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1730 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1732 GetDigitsArray()->GetTree()->SetName(treeName);
1733 GetDigitsArray()->GetTree()->AutoSave();
1740 //__________________________________________________________________
1741 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1744 //-----------------------------------------------------------
1745 // summable digits - 16 bit "ADC", no noise, no saturation
1746 //-----------------------------------------------------------
1748 //----------------------------------------------------
1749 // Loop over all sectors for a single event
1750 //----------------------------------------------------
1753 if(fDefaults == 0) SetDefaults();
1754 GenerNoise(500000); //create table with noise
1755 //setup TPCDigitsArray
1757 if(GetDigitsArray()) delete GetDigitsArray();
1759 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1760 arr->SetClass("AliSimDigits");
1761 arr->Setup(fTPCParam);
1762 // Note that methods arr->MakeTree have different signatures
1763 if (gAlice->GetTreeSFile()) {
1764 arr->MakeTree(gAlice->GetTreeSFile());
1766 arr->MakeTree(fDigitsFile);
1768 SetDigitsArray(arr);
1770 cerr<<"Digitizing TPC -- summable digits...\n";
1772 fDigitsSwitch=1; // summable digits
1774 // set zero suppression to "0"
1776 fTPCParam->SetZeroSup(0);
1778 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1785 sprintf(treeName,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1787 GetDigitsArray()->GetTree()->SetName(treeName);
1788 GetDigitsArray()->GetTree()->AutoSave();
1794 //__________________________________________________________________
1795 void AliTPC::Hits2SDigits()
1798 //-----------------------------------------------------------
1799 // summable digits - 16 bit "ADC", no noise, no saturation
1800 //-----------------------------------------------------------
1802 //----------------------------------------------------
1803 // Loop over all sectors for a single event
1804 //----------------------------------------------------
1805 //MI change - for pp run
1806 Int_t eventnumber = gAlice->GetEvNumber();
1808 if(fDefaults == 0) SetDefaults();
1809 GenerNoise(500000); //create table with noise
1811 //setup TPCDigitsArray
1813 if(GetDigitsArray()) delete GetDigitsArray();
1815 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1816 arr->SetClass("AliSimDigits");
1817 arr->Setup(fTPCParam);
1818 // Note that methods arr->MakeTree have different signatures
1819 if (gAlice->GetTreeSFile()) {
1820 arr->MakeTree(gAlice->GetTreeSFile());
1822 arr->MakeTree(fDigitsFile);
1824 SetDigitsArray(arr);
1826 cerr<<"Digitizing TPC -- summable digits...\n";
1828 // fDigitsSwitch=1; // summable digits -for the moment direct
1830 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1836 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1838 GetDigitsArray()->GetTree()->SetName(treeName);
1839 GetDigitsArray()->GetTree()->AutoSave();
1844 //_____________________________________________________________________________
1845 void AliTPC::Hits2DigitsSector(Int_t isec)
1847 //-------------------------------------------------------------------
1848 // TPC conversion from hits to digits.
1849 //-------------------------------------------------------------------
1851 //-----------------------------------------------------------------
1852 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1853 //-----------------------------------------------------------------
1855 //-------------------------------------------------------
1856 // Get the access to the track hits
1857 //-------------------------------------------------------
1859 // check if the parameters are set - important if one calls this method
1860 // directly, not from the Hits2Digits
1862 if(fDefaults == 0) SetDefaults();
1864 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1865 Stat_t ntracks = tH->GetEntries();
1869 //-------------------------------------------
1870 // Only if there are any tracks...
1871 //-------------------------------------------
1875 //printf("*** Processing sector number %d ***\n",isec);
1877 Int_t nrows =fTPCParam->GetNRow(isec);
1879 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1881 MakeSector(isec,nrows,tH,ntracks,row);
1883 //--------------------------------------------------------
1884 // Digitize this sector, row by row
1885 // row[i] is the pointer to the TObjArray of AliTPCFastVectors,
1886 // each one containing electrons accepted on this
1887 // row, assigned into tracks
1888 //--------------------------------------------------------
1892 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1894 for (i=0;i<nrows;i++){
1896 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1898 DigitizeRow(i,isec,row);
1900 fDigitsArray->StoreRow(isec,i);
1902 Int_t ndig = dig->GetDigitSize();
1903 if (gDebug > 10) printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1905 fDigitsArray->ClearRow(isec,i);
1908 } // end of the sector digitization
1910 for(i=0;i<nrows+2;i++){
1915 delete [] row; // delete the array of pointers to TObjArray-s
1919 } // end of Hits2DigitsSector
1922 //_____________________________________________________________________________
1923 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1925 //-----------------------------------------------------------
1926 // Single row digitization, coupling from the neighbouring
1927 // rows taken into account
1928 //-----------------------------------------------------------
1930 //-----------------------------------------------------------------
1931 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1932 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1933 //-----------------------------------------------------------------
1936 Float_t zerosup = fTPCParam->GetZeroSup();
1937 // Int_t nrows =fTPCParam->GetNRow(isec);
1938 fCurrentIndex[1]= isec;
1941 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1942 Int_t nofTbins = fTPCParam->GetMaxTBin();
1943 Int_t indexRange[4];
1945 // Integrated signal for this row
1946 // and a single track signal
1949 AliTPCFastMatrix *m1 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // integrated
1950 AliTPCFastMatrix *m2 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // single
1952 AliTPCFastMatrix &total = *m1;
1954 // Array of pointers to the label-signal list
1956 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1957 Float_t **pList = new Float_t* [nofDigits];
1961 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1965 //Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1966 //Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1969 for (Int_t row= row1;row<=row2;row++){
1970 Int_t nTracks= rows[row]->GetEntries();
1971 for (i1=0;i1<nTracks;i1++){
1972 fCurrentIndex[2]= row;
1973 fCurrentIndex[3]=irow+1;
1975 m2->Zero(); // clear single track signal matrix
1976 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1977 GetList(trackLabel,nofPads,m2,indexRange,pList);
1979 else GetSignal(rows[row],i1,0,m1,indexRange);
1985 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1987 Float_t fzerosup = zerosup+0.5;
1988 for(Int_t it=0;it<nofTbins;it++){
1989 Float_t *pq = &(total.UncheckedAt(0,it));
1990 for(Int_t ip=0;ip<nofPads;ip++){
1994 if(fDigitsSwitch == 0){
1996 if(q <=fzerosup) continue; // do not fill zeros
1998 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
2003 if(q <= 0.) continue; // do not fill zeros
2004 if(q>2000.) q=2000.;
2010 // "real" signal or electronic noise (list = -1)?
2013 for(Int_t j1=0;j1<3;j1++){
2014 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
2019 <A NAME="AliDigits"></A>
2020 using of AliDigits object
2023 dig->SetDigitFast((Short_t)q,it,ip);
2024 if (fDigitsArray->IsSimulated())
2026 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
2027 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
2028 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
2032 } // end of loop over time buckets
2033 } // end of lop over pads
2036 // This row has been digitized, delete nonused stuff
2039 for(lp=0;lp<nofDigits;lp++){
2040 if(pList[lp]) delete [] pList[lp];
2049 } // end of DigitizeRow
2051 //_____________________________________________________________________________
2053 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
2054 AliTPCFastMatrix *m1, AliTPCFastMatrix *m2,Int_t *indexRange)
2057 //---------------------------------------------------------------
2058 // Calculates 2-D signal (pad,time) for a single track,
2059 // returns a pointer to the signal matrix and the track label
2060 // No digitization is performed at this level!!!
2061 //---------------------------------------------------------------
2063 //-----------------------------------------------------------------
2064 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2065 // Modified: Marian Ivanov
2066 //-----------------------------------------------------------------
2068 AliTPCFastVector *tv;
2070 tv = (AliTPCFastVector*)p1->At(ntr); // pointer to a track
2071 AliTPCFastVector &v = *tv;
2073 Float_t label = v(0);
2074 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1)-1)/2;
2076 Int_t nElectrons = (tv->GetNrows()-1)/4;
2077 indexRange[0]=9999; // min pad
2078 indexRange[1]=-1; // max pad
2079 indexRange[2]=9999; //min time
2080 indexRange[3]=-1; // max time
2082 AliTPCFastMatrix &signal = *m1;
2083 AliTPCFastMatrix &total = *m2;
2085 // Loop over all electrons
2087 for(Int_t nel=0; nel<nElectrons; nel++){
2089 Float_t aval = v(idx+4);
2090 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
2091 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
2092 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
2094 Int_t *index = fTPCParam->GetResBin(0);
2095 Float_t *weight = & (fTPCParam->GetResWeight(0));
2097 if (n>0) for (Int_t i =0; i<n; i++){
2098 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
2101 Int_t time=index[2];
2102 Float_t qweight = *(weight)*eltoadcfac;
2104 if (m1!=0) signal.UncheckedAt(pad,time)+=qweight;
2105 total.UncheckedAt(pad,time)+=qweight;
2106 if (indexRange[0]>pad) indexRange[0]=pad;
2107 if (indexRange[1]<pad) indexRange[1]=pad;
2108 if (indexRange[2]>time) indexRange[2]=time;
2109 if (indexRange[3]<time) indexRange[3]=time;
2116 } // end of loop over electrons
2118 return label; // returns track label when finished
2121 //_____________________________________________________________________________
2122 void AliTPC::GetList(Float_t label,Int_t np,AliTPCFastMatrix *m,
2123 Int_t *indexRange, Float_t **pList)
2125 //----------------------------------------------------------------------
2126 // Updates the list of tracks contributing to digits for a given row
2127 //----------------------------------------------------------------------
2129 //-----------------------------------------------------------------
2130 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2131 //-----------------------------------------------------------------
2133 AliTPCFastMatrix &signal = *m;
2135 // lop over nonzero digits
2137 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
2138 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
2141 // accept only the contribution larger than 500 electrons (1/2 s_noise)
2143 if(signal(ip,it)<0.5) continue;
2146 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
2148 if(!pList[globalIndex]){
2151 // Create new list (6 elements - 3 signals and 3 labels),
2154 pList[globalIndex] = new Float_t [6];
2158 *pList[globalIndex] = -1.;
2159 *(pList[globalIndex]+1) = -1.;
2160 *(pList[globalIndex]+2) = -1.;
2161 *(pList[globalIndex]+3) = -1.;
2162 *(pList[globalIndex]+4) = -1.;
2163 *(pList[globalIndex]+5) = -1.;
2166 *pList[globalIndex] = label;
2167 *(pList[globalIndex]+3) = signal(ip,it);
2171 // check the signal magnitude
2173 Float_t highest = *(pList[globalIndex]+3);
2174 Float_t middle = *(pList[globalIndex]+4);
2175 Float_t lowest = *(pList[globalIndex]+5);
2178 // compare the new signal with already existing list
2181 if(signal(ip,it)<lowest) continue; // neglect this track
2185 if (signal(ip,it)>highest){
2186 *(pList[globalIndex]+5) = middle;
2187 *(pList[globalIndex]+4) = highest;
2188 *(pList[globalIndex]+3) = signal(ip,it);
2190 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2191 *(pList[globalIndex]+1) = *pList[globalIndex];
2192 *pList[globalIndex] = label;
2194 else if (signal(ip,it)>middle){
2195 *(pList[globalIndex]+5) = middle;
2196 *(pList[globalIndex]+4) = signal(ip,it);
2198 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2199 *(pList[globalIndex]+1) = label;
2202 *(pList[globalIndex]+5) = signal(ip,it);
2203 *(pList[globalIndex]+2) = label;
2207 } // end of loop over pads
2208 } // end of loop over time bins
2213 //___________________________________________________________________
2214 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
2215 Stat_t ntracks,TObjArray **row)
2218 //-----------------------------------------------------------------
2219 // Prepares the sector digitization, creates the vectors of
2220 // tracks for each row of this sector. The track vector
2221 // contains the track label and the position of electrons.
2222 //-----------------------------------------------------------------
2224 //-----------------------------------------------------------------
2225 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2226 //-----------------------------------------------------------------
2228 Float_t gasgain = fTPCParam->GetGasGain();
2232 AliTPChit *tpcHit; // pointer to a sigle TPC hit
2235 if (fHitType>1) branch = TH->GetBranch("TPC2");
2236 else branch = TH->GetBranch("TPC");
2239 //----------------------------------------------
2240 // Create TObjArray-s, one for each row,
2241 // each TObjArray will store the AliTPCFastVectors
2242 // of electrons, one AliTPCFastVectors per each track.
2243 //----------------------------------------------
2245 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
2246 AliTPCFastVector **tracks = new AliTPCFastVector* [nrows+2]; //pointers to the track vectors
2248 for(i=0; i<nrows+2; i++){
2249 row[i] = new TObjArray;
2256 //--------------------------------------------------------------------
2257 // Loop over tracks, the "track" contains the full history
2258 //--------------------------------------------------------------------
2260 Int_t previousTrack,currentTrack;
2261 previousTrack = -1; // nothing to store so far!
2263 for(Int_t track=0;track<ntracks;track++){
2264 Bool_t isInSector=kTRUE;
2266 isInSector = TrackInVolume(isec,track);
2267 if (!isInSector) continue;
2269 branch->GetEntry(track); // get next track
2273 tpcHit = (AliTPChit*)FirstHit(-1);
2275 //--------------------------------------------------------------
2277 //--------------------------------------------------------------
2282 Int_t sector=tpcHit->fSector; // sector number
2284 tpcHit = (AliTPChit*) NextHit();
2288 currentTrack = tpcHit->Track(); // track number
2291 if(currentTrack != previousTrack){
2293 // store already filled fTrack
2295 for(i=0;i<nrows+2;i++){
2296 if(previousTrack != -1){
2297 if(nofElectrons[i]>0){
2298 AliTPCFastVector &v = *tracks[i];
2299 v(0) = previousTrack;
2300 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2301 row[i]->Add(tracks[i]);
2304 delete tracks[i]; // delete empty AliTPCFastVector
2310 tracks[i] = new AliTPCFastVector(481); // AliTPCFastVectors for the next fTrack
2312 } // end of loop over rows
2314 previousTrack=currentTrack; // update track label
2317 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2319 //---------------------------------------------------
2320 // Calculate the electron attachment probability
2321 //---------------------------------------------------
2324 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
2325 /fTPCParam->GetDriftV();
2327 Float_t attProb = fTPCParam->GetAttCoef()*
2328 fTPCParam->GetOxyCont()*time; // fraction!
2330 //-----------------------------------------------
2331 // Loop over electrons
2332 //-----------------------------------------------
2335 for(Int_t nel=0;nel<qI;nel++){
2336 // skip if electron lost due to the attachment
2337 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2342 // protection for the nonphysical avalanche size (10**6 maximum)
2344 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2345 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2348 TransportElectron(xyz,index);
2350 fTPCParam->GetPadRow(xyz,index);
2351 // row 0 - cross talk from the innermost row
2352 // row fNRow+1 cross talk from the outermost row
2353 rowNumber = index[2]+1;
2354 //transform position to local digit coordinates
2355 //relative to nearest pad row
2356 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2358 if (isec <fTPCParam->GetNInnerSector()) {
2359 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2360 y1 = fTPCParam->GetYInner(rowNumber);
2363 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2364 y1 = fTPCParam->GetYOuter(rowNumber);
2366 // gain inefficiency at the wires edges - linear
2369 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));
2371 nofElectrons[rowNumber]++;
2372 //----------------------------------
2373 // Expand vector if necessary
2374 //----------------------------------
2375 if(nofElectrons[rowNumber]>120){
2376 Int_t range = tracks[rowNumber]->GetNrows();
2377 if((nofElectrons[rowNumber])>(range-1)/4){
2379 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
2383 AliTPCFastVector &v = *tracks[rowNumber];
2384 Int_t idx = 4*nofElectrons[rowNumber]-3;
2385 Real_t * position = &(((AliTPCFastVector&)v).UncheckedAt(idx)); //make code faster
2386 memcpy(position,xyz,4*sizeof(Float_t));
2388 } // end of loop over electrons
2390 tpcHit = (AliTPChit*)NextHit();
2392 } // end of loop over hits
2393 } // end of loop over tracks
2396 // store remaining track (the last one) if not empty
2399 for(i=0;i<nrows+2;i++){
2400 if(nofElectrons[i]>0){
2401 AliTPCFastVector &v = *tracks[i];
2402 v(0) = previousTrack;
2403 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2404 row[i]->Add(tracks[i]);
2413 delete [] nofElectrons;
2416 } // end of MakeSector
2419 //_____________________________________________________________________________
2423 // Initialise TPC detector after definition of geometry
2428 printf("\n%s: ",ClassName());
2429 for(i=0;i<35;i++) printf("*");
2430 printf(" TPC_INIT ");
2431 for(i=0;i<35;i++) printf("*");
2432 printf("\n%s: ",ClassName());
2434 for(i=0;i<80;i++) printf("*");
2439 //_____________________________________________________________________________
2440 void AliTPC::MakeBranch(Option_t* option, const char *file)
2443 // Create Tree branches for the TPC.
2445 Int_t buffersize = 4000;
2446 char branchname[10];
2447 sprintf(branchname,"%s",GetName());
2449 AliDetector::MakeBranch(option,file);
2451 const char *d = strstr(option,"D");
2453 if (fDigits && gAlice->TreeD() && d) {
2454 MakeBranchInTree(gAlice->TreeD(),
2455 branchname, &fDigits, buffersize, file);
2458 if (fHitType>1) MakeBranch2(option,file); // MI change 14.09.2000
2461 //_____________________________________________________________________________
2462 void AliTPC::ResetDigits()
2465 // Reset number of digits and the digits array for this detector
2468 if (fDigits) fDigits->Clear();
2471 //_____________________________________________________________________________
2472 void AliTPC::SetSecAL(Int_t sec)
2474 //---------------------------------------------------
2475 // Activate/deactivate selection for lower sectors
2476 //---------------------------------------------------
2478 //-----------------------------------------------------------------
2479 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2480 //-----------------------------------------------------------------
2485 //_____________________________________________________________________________
2486 void AliTPC::SetSecAU(Int_t sec)
2488 //----------------------------------------------------
2489 // Activate/deactivate selection for upper sectors
2490 //---------------------------------------------------
2492 //-----------------------------------------------------------------
2493 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2494 //-----------------------------------------------------------------
2499 //_____________________________________________________________________________
2500 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2502 //----------------------------------------
2503 // Select active lower sectors
2504 //----------------------------------------
2506 //-----------------------------------------------------------------
2507 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2508 //-----------------------------------------------------------------
2518 //_____________________________________________________________________________
2519 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2520 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2521 Int_t s11 , Int_t s12)
2523 //--------------------------------
2524 // Select active upper sectors
2525 //--------------------------------
2527 //-----------------------------------------------------------------
2528 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2529 //-----------------------------------------------------------------
2545 //_____________________________________________________________________________
2546 void AliTPC::SetSens(Int_t sens)
2549 //-------------------------------------------------------------
2550 // Activates/deactivates the sensitive strips at the center of
2551 // the pad row -- this is for the space-point resolution calculations
2552 //-------------------------------------------------------------
2554 //-----------------------------------------------------------------
2555 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2556 //-----------------------------------------------------------------
2562 void AliTPC::SetSide(Float_t side=0.)
2564 // choice of the TPC side
2569 //____________________________________________________________________________
2570 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2571 Float_t p2,Float_t p3)
2574 // gax mixture definition
2588 //_____________________________________________________________________________
2590 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2593 // electron transport taking into account:
2595 // 2.ExB at the wires
2596 // 3. nonisochronity
2598 // xyz and index must be already transformed to system 1
2601 fTPCParam->Transform1to2(xyz,index);
2604 Float_t driftl=xyz[2];
2605 if(driftl<0.01) driftl=0.01;
2606 driftl=TMath::Sqrt(driftl);
2607 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2608 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2609 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2610 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2611 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2615 if (fTPCParam->GetMWPCReadout()==kTRUE){
2616 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2617 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2619 //add nonisochronity (not implemented yet)
2622 ClassImp(AliTPCdigit)
2624 //_____________________________________________________________________________
2625 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2629 // Creates a TPC digit object
2631 fSector = digits[0];
2632 fPadRow = digits[1];
2635 fSignal = digits[4];
2641 //_____________________________________________________________________________
2642 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2646 // Creates a TPC hit object
2657 //________________________________________________________________________
2658 // Additional code because of the AliTPCTrackHitsV2
2660 void AliTPC::MakeBranch2(Option_t *option,const char *file)
2663 // Create a new branch in the current Root Tree
2664 // The branch of fHits is automatically split
2665 // MI change 14.09.2000
2666 if (fHitType<2) return;
2667 char branchname[10];
2668 sprintf(branchname,"%s2",GetName());
2670 // Get the pointer to the header
2671 const char *cH = strstr(option,"H");
2673 if (fTrackHits && gAlice->TreeH() && cH && fHitType&4) {
2674 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2675 // gAlice->TreeH(),fBufferSize,99);
2676 //TBranch * branch =
2677 gAlice->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2680 // gAlice->TreeH()->GetListOfBranches()->Add(branch);
2682 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2683 const char folder [] = "RunMC/Event/Data";
2685 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2686 Publish(folder,&fTrackHits,branchname);
2688 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2689 TDirectory *wd = gDirectory;
2691 TIter next( b->GetListOfBranches());
2692 while ((b=(TBranch*)next())) {
2697 cout << "Diverting branch " << branchname << " to file " << file << endl;
2701 if (fTrackHitsOld && gAlice->TreeH() && cH && fHitType&2) {
2702 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2703 gAlice->TreeH(),fBufferSize,99);
2704 //TBranch * branch = gAlice->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2707 gAlice->TreeH()->GetListOfBranches()->Add(branch);
2709 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2710 const char folder [] = "RunMC/Event/Data";
2712 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2713 Publish(folder,&fTrackHitsOld,branchname);
2715 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2716 TDirectory *wd = gDirectory;
2718 TIter next( b->GetListOfBranches());
2719 while ((b=(TBranch*)next())) {
2724 cout << "Diverting branch " << branchname << " to file " << file << endl;
2729 void AliTPC::SetTreeAddress()
2731 if (fHitType<=1) AliDetector::SetTreeAddress();
2732 if (fHitType>1) SetTreeAddress2();
2735 void AliTPC::SetTreeAddress2()
2738 // Set branch address for the TrackHits Tree
2741 char branchname[20];
2742 sprintf(branchname,"%s2",GetName());
2744 // Branch address for hit tree
2745 TTree *treeH = gAlice->TreeH();
2746 if ((treeH)&&(fHitType&4)) {
2747 branch = treeH->GetBranch(branchname);
2748 if (branch) branch->SetAddress(&fTrackHits);
2750 if ((treeH)&&(fHitType&2)) {
2751 branch = treeH->GetBranch(branchname);
2752 if (branch) branch->SetAddress(&fTrackHitsOld);
2754 //set address to TREETR
2755 TTree *treeTR = gAlice->TreeTR();
2756 if (treeTR && fTrackReferences) {
2757 branch = treeTR->GetBranch(GetName());
2758 if (branch) branch->SetAddress(&fTrackReferences);
2764 void AliTPC::FinishPrimary()
2766 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2767 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2771 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2774 // add hit to the list
2777 int primary = gAlice->GetPrimary(track);
2778 gAlice->Particle(primary)->SetBit(kKeepBit);
2782 gAlice->FlagTrack(track);
2784 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2785 //if (hit->fTrack!=rtrack)
2786 // cout<<"bad track number\n";
2787 if (fTrackHits && fHitType&4)
2788 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2789 hits[1],hits[2],(Int_t)hits[3]);
2790 if (fTrackHitsOld &&fHitType&2 )
2791 fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2792 hits[1],hits[2],(Int_t)hits[3]);
2796 void AliTPC::ResetHits()
2798 if (fHitType&1) AliDetector::ResetHits();
2799 if (fHitType>1) ResetHits2();
2802 void AliTPC::ResetHits2()
2806 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2807 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2811 AliHit* AliTPC::FirstHit(Int_t track)
2813 if (fHitType>1) return FirstHit2(track);
2814 return AliDetector::FirstHit(track);
2816 AliHit* AliTPC::NextHit()
2818 if (fHitType>1) return NextHit2();
2820 return AliDetector::NextHit();
2823 AliHit* AliTPC::FirstHit2(Int_t track)
2826 // Initialise the hit iterator
2827 // Return the address of the first hit for track
2828 // If track>=0 the track is read from disk
2829 // while if track<0 the first hit of the current
2830 // track is returned
2833 gAlice->ResetHits();
2834 gAlice->TreeH()->GetEvent(track);
2837 if (fTrackHits && fHitType&4) {
2838 fTrackHits->First();
2839 return fTrackHits->GetHit();
2841 if (fTrackHitsOld && fHitType&2) {
2842 fTrackHitsOld->First();
2843 return fTrackHitsOld->GetHit();
2849 AliHit* AliTPC::NextHit2()
2852 //Return the next hit for the current track
2855 if (fTrackHitsOld && fHitType&2) {
2856 fTrackHitsOld->Next();
2857 return fTrackHitsOld->GetHit();
2861 return fTrackHits->GetHit();
2867 void AliTPC::LoadPoints(Int_t)
2871 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2874 if(fHitType==1) AliDetector::LoadPoints(a);
2875 else LoadPoints2(a);
2882 void AliTPC::RemapTrackHitIDs(Int_t *map)
2884 if (!fTrackHits) return;
2886 if (fTrackHitsOld && fHitType&2){
2887 AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2888 for (UInt_t i=0;i<arr->GetSize();i++){
2889 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2890 info->fTrackID = map[info->fTrackID];
2893 if (fTrackHitsOld && fHitType&4){
2894 TClonesArray * arr = fTrackHits->GetArray();;
2895 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2896 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2897 info->fTrackID = map[info->fTrackID];
2902 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2904 //return bool information - is track in given volume
2905 //load only part of the track information
2906 //return true if current track is in volume
2909 if (fTrackHitsOld && fHitType&2) {
2910 TBranch * br = gAlice->TreeH()->GetBranch("fTrackHitsInfo");
2911 br->GetEvent(track);
2912 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2913 for (UInt_t j=0;j<ar->GetSize();j++){
2914 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2918 if (fTrackHits && fHitType&4) {
2919 TBranch * br1 = gAlice->TreeH()->GetBranch("fVolumes");
2920 TBranch * br2 = gAlice->TreeH()->GetBranch("fNVolumes");
2921 br2->GetEvent(track);
2922 br1->GetEvent(track);
2923 Int_t *volumes = fTrackHits->GetVolumes();
2924 Int_t nvolumes = fTrackHits->GetNVolumes();
2925 if (!volumes && nvolumes>0) {
2926 printf("Problematic track\t%d\t%d",track,nvolumes);
2929 for (Int_t j=0;j<nvolumes; j++)
2930 if (volumes[j]==id) return kTRUE;
2934 TBranch * br = gAlice->TreeH()->GetBranch("fSector");
2935 br->GetEvent(track);
2936 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2937 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2944 //_____________________________________________________________________________
2945 void AliTPC::LoadPoints2(Int_t)
2948 // Store x, y, z of all hits in memory
2950 if (fTrackHits == 0 && fTrackHitsOld==0) return;
2953 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2954 if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2956 if (nhits == 0) return;
2957 Int_t tracks = gAlice->GetNtrack();
2958 if (fPoints == 0) fPoints = new TObjArray(tracks);
2961 Int_t *ntrk=new Int_t[tracks];
2962 Int_t *limi=new Int_t[tracks];
2963 Float_t **coor=new Float_t*[tracks];
2964 for(Int_t i=0;i<tracks;i++) {
2970 AliPoints *points = 0;
2973 Int_t chunk=nhits/4+1;
2975 // Loop over all the hits and store their position
2977 ahit = FirstHit2(-1);
2979 trk=ahit->GetTrack();
2980 if(ntrk[trk]==limi[trk]) {
2982 // Initialise a new track
2983 fp=new Float_t[3*(limi[trk]+chunk)];
2985 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2986 delete [] coor[trk];
2993 fp[3*ntrk[trk] ] = ahit->X();
2994 fp[3*ntrk[trk]+1] = ahit->Y();
2995 fp[3*ntrk[trk]+2] = ahit->Z();
3003 for(trk=0; trk<tracks; ++trk) {
3005 points = new AliPoints();
3006 points->SetMarkerColor(GetMarkerColor());
3007 points->SetMarkerSize(GetMarkerSize());
3008 points->SetDetector(this);
3009 points->SetParticle(trk);
3010 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
3011 fPoints->AddAt(points,trk);
3012 delete [] coor[trk];
3022 //_____________________________________________________________________________
3023 void AliTPC::LoadPoints3(Int_t)
3026 // Store x, y, z of all hits in memory
3027 // - only intersection point with pad row
3028 if (fTrackHits == 0) return;
3030 Int_t nhits = fTrackHits->GetEntriesFast();
3031 if (nhits == 0) return;
3032 Int_t tracks = gAlice->GetNtrack();
3033 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
3034 fPoints->Expand(2*tracks);
3037 Int_t *ntrk=new Int_t[tracks];
3038 Int_t *limi=new Int_t[tracks];
3039 Float_t **coor=new Float_t*[tracks];
3040 for(Int_t i=0;i<tracks;i++) {
3046 AliPoints *points = 0;
3049 Int_t chunk=nhits/4+1;
3051 // Loop over all the hits and store their position
3053 ahit = FirstHit2(-1);
3054 //for (Int_t hit=0;hit<nhits;hit++) {
3058 // ahit = (AliHit*)fHits->UncheckedAt(hit);
3059 trk=ahit->GetTrack();
3060 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
3061 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
3062 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
3063 if (currentrow!=lastrow){
3064 lastrow = currentrow;
3065 //later calculate intersection point
3066 if(ntrk[trk]==limi[trk]) {
3068 // Initialise a new track
3069 fp=new Float_t[3*(limi[trk]+chunk)];
3071 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
3072 delete [] coor[trk];
3079 fp[3*ntrk[trk] ] = ahit->X();
3080 fp[3*ntrk[trk]+1] = ahit->Y();
3081 fp[3*ntrk[trk]+2] = ahit->Z();
3088 for(trk=0; trk<tracks; ++trk) {
3090 points = new AliPoints();
3091 points->SetMarkerColor(GetMarkerColor()+1);
3092 points->SetMarkerStyle(5);
3093 points->SetMarkerSize(0.2);
3094 points->SetDetector(this);
3095 points->SetParticle(trk);
3096 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
3097 points->SetPolyMarker(ntrk[trk],coor[trk],30);
3098 fPoints->AddAt(points,tracks+trk);
3099 delete [] coor[trk];
3110 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
3114 //fill clones array with intersection of current point with the
3119 const Int_t kcmaxhits=30000;
3120 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
3121 AliTPCFastVector & xxx = *xxxx;
3122 Int_t maxhits = kcmaxhits;
3125 AliTPChit * tpcHit=0;
3126 tpcHit = (AliTPChit*)FirstHit2(-1);
3127 Int_t currentIndex=0;
3128 Int_t lastrow=-1; //last writen row
3131 if (tpcHit==0) continue;
3132 sector=tpcHit->fSector; // sector number
3133 ipart=tpcHit->Track();
3137 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
3138 Int_t index[3]={1,sector,0};
3139 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
3140 if (currentrow<0) continue;
3141 if (lastrow<0) lastrow=currentrow;
3142 if (currentrow==lastrow){
3143 if ( currentIndex>=maxhits){
3145 xxx.ResizeTo(4*maxhits);
3147 xxx(currentIndex*4)=x[0];
3148 xxx(currentIndex*4+1)=x[1];
3149 xxx(currentIndex*4+2)=x[2];
3150 xxx(currentIndex*4+3)=tpcHit->fQ;
3154 if (currentIndex>2){
3166 for (Int_t index=0;index<currentIndex;index++){
3168 x=x2=x3=x4=xxx(index*4);
3176 sumy+=xxx(index*4+1);
3177 sumxy+=xxx(index*4+1)*x;
3178 sumx2y+=xxx(index*4+1)*x2;
3179 sumz+=xxx(index*4+2);
3180 sumxz+=xxx(index*4+2)*x;
3181 sumx2z+=xxx(index*4+2)*x2;
3182 sumq+=xxx(index*4+3);
3184 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
3185 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
3186 sumx2*(sumx*sumx3-sumx2*sumx2);
3188 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
3189 sumx2*(sumxy*sumx3-sumx2y*sumx2);
3190 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
3191 sumx2*(sumxz*sumx3-sumx2z*sumx2);
3193 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
3194 sumx2*(sumx*sumx2y-sumx2*sumxy);
3195 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
3196 sumx2*(sumx*sumx2z-sumx2*sumxz);
3198 Float_t y=detay/det+centralPad;
3199 Float_t z=detaz/det;
3200 Float_t by=detby/det; //y angle
3201 Float_t bz=detbz/det; //z angle
3202 sumy/=Float_t(currentIndex);
3203 sumz/=Float_t(currentIndex);
3205 AliComplexCluster cl;
3211 cl.fTracks[0]=ipart;
3213 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
3214 if (row!=0) row->InsertCluster(&cl);
3217 } //end of calculating cluster for given row
3219 } // end of loop over hits
3223 //_______________________________________________________________________________
3224 void AliTPC::Digits2Reco(Int_t firstevent,Int_t lastevent)
3226 // produces rec points from digits and writes them on the root file
3227 // AliTPCclusters.root
3229 TDirectory *cwd = gDirectory;
3232 AliTPCParamSR *dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60");
3234 printf("You are running 2 pad-length geom hits with 3 pad-length geom digits\n");
3236 dig = new AliTPCParamSR();
3240 dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");
3243 printf("No TPC parameters found\n");
3248 cout<<"AliTPC::Digits2Reco: TPC parameteres have been set"<<endl;
3250 if(!gSystem->Getenv("CONFIG_FILE")){
3251 out=TFile::Open("AliTPCclusters.root","recreate");
3257 tmp1=gSystem->Getenv("CONFIG_FILE_PREFIX");
3258 tmp2=gSystem->Getenv("CONFIG_OUTDIR");
3259 sprintf(tmp3,"%s%s/AliTPCclusters.root",tmp1,tmp2);
3260 out=TFile::Open(tmp3,"recreate");
3264 cout<<"AliTPC::Digits2Reco - determination of rec points begins"<<endl;
3267 for(Int_t iev=firstevent;iev<lastevent+1;iev++){
3269 printf("Processing event %d\n",iev);
3270 Digits2Clusters(out,iev);
3272 cout<<"AliTPC::Digits2Reco - determination of rec points ended"<<endl;
3281 ////////////////////////////////////////////////////////////////////////
3282 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
3284 // load TPC paarmeters from a given file or create new if the object
3285 // is not found there
3288 sprintf(paramName,"75x40_100x60_150x60");
3289 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
3291 cout<<"TPC parameters "<<paramName<<" found."<<endl;
3293 cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3295 paramTPC = new AliTPCParamSR;
3299 // the older version of parameters can be accessed with this code.
3300 // In some cases, we have old parameters saved in the file but
3301 // digits were created with new parameters, it can be distinguish
3302 // by the name of TPC TreeD. The code here is just for the case
3303 // we would need to compare with old data, uncomment it if needed.
3305 // char paramName[50];
3306 // sprintf(paramName,"75x40_100x60");
3307 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
3309 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
3311 // sprintf(paramName,"75x40_100x60_150x60");
3312 // paramTPC=(AliTPCParam*)in->Get(paramName);
3314 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
3316 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3318 // paramTPC = new AliTPCParamSR;