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.58 2002/05/27 14:33:14 hristov
19 The new class AliTrackReference used (M.Ivanov)
21 Revision 1.57 2002/05/07 17:23:11 kowal2
22 Linear gain inefficiency instead of the step one at the wire edges.
24 Revision 1.56 2002/04/04 16:26:33 kowal2
25 Digits (Sdigits) go to separate files now.
27 Revision 1.55 2002/03/29 06:57:45 kowal2
28 Restored backward compatibility to use the hits from Dec. 2000 production.
30 Revision 1.54 2002/03/18 17:59:13 kowal2
31 Chnges in the pad geometry - 3 pad lengths introduced.
33 Revision 1.53 2002/02/25 11:02:56 kowal2
34 Changes towards speeding up the code. Thanks to Marian Ivanov.
36 Revision 1.52 2002/02/18 09:26:09 kowal2
37 Removed compiler warning
39 Revision 1.51 2002/01/21 17:13:21 kowal2
40 New track hits using root containers. Setting active sectors added.
42 Revision 1.50 2001/12/06 14:16:19 kowal2
45 Revision 1.49 2001/11/30 11:55:37 hristov
46 Noise table created in Hits2SDigits (M.Ivanov)
48 Revision 1.48 2001/11/24 16:10:21 kowal2
51 Revision 1.47 2001/11/19 10:25:34 kowal2
52 Nearest integer instead of integer when converting to ADC counts
54 Revision 1.46 2001/11/07 06:47:12 kowal2
57 Revision 1.45 2001/11/03 13:33:48 kowal2
58 Updated algorithms in Hits2SDigits, SDigits2Digits,
62 Revision 1.44 2001/08/30 09:28:48 hristov
63 TTree names are explicitly set via SetName(name) and then Write() is called
65 Revision 1.43 2001/07/28 12:02:54 hristov
66 Branch split level set to 99
68 Revision 1.42 2001/07/28 11:38:52 hristov
69 Loop variable declared once
71 Revision 1.41 2001/07/28 10:53:50 hristov
72 Digitisation done according to the general scheme (M.Ivanov)
74 Revision 1.40 2001/07/27 13:03:14 hristov
75 Default Branch split level set to 99
77 Revision 1.39 2001/07/26 09:09:34 kowal2
78 Hits2Reco method added
80 Revision 1.38 2001/07/20 14:32:43 kowal2
81 Processing of many events possible now
83 Revision 1.37 2001/06/12 07:17:18 kowal2
84 Hits2SDigits method implemented (summable digits)
86 Revision 1.36 2001/05/16 14:57:25 alibrary
87 New files for folders and Stack
89 Revision 1.35 2001/05/08 16:02:22 kowal2
90 Updated material specifications
92 Revision 1.34 2001/05/08 15:00:15 hristov
93 Corrections for tracking in arbitrary magnenetic field. Changes towards a concept of global Alice track. Back propagation of reconstructed tracks (Yu.Belikov)
95 Revision 1.33 2001/04/03 12:40:43 kowal2
98 Revision 1.32 2001/03/12 17:47:36 hristov
99 Changes needed on Sun with CC 5.0
101 Revision 1.31 2001/03/12 08:21:50 kowal2
102 Corrected C++ bug in the material definitions
104 Revision 1.30 2001/03/01 17:34:47 kowal2
105 Correction due to the accuracy problem
107 Revision 1.29 2001/02/28 16:34:40 kowal2
108 Protection against nonphysical values of the avalanche size,
111 Revision 1.28 2001/01/26 19:57:19 hristov
112 Major upgrade of AliRoot code
114 Revision 1.27 2001/01/13 17:29:33 kowal2
115 Sun compiler correction
117 Revision 1.26 2001/01/10 07:59:43 kowal2
118 Corrections to load points from the noncompressed hits.
120 Revision 1.25 2000/11/02 07:25:31 kowal2
121 Changes due to the new hit structure.
124 Revision 1.24 2000/10/05 16:06:09 kowal2
125 Forward declarations. Changes due to a new class AliComplexCluster.
127 Revision 1.23 2000/10/02 21:28:18 fca
128 Removal of useless dependecies via forward declarations
130 Revision 1.22 2000/07/10 20:57:39 hristov
131 Update of TPC code and macros by M.Kowalski
133 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
134 Changes to obey the coding rules
136 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
137 Splitted from AliTPCtracking
139 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
140 Changed parameter settings
142 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
144 Defaults loaded automatically (hard-wired)
145 Optional parameters can be set via macro called in the constructor
147 Revision 1.19 2000/04/18 19:00:59 fca
148 Small bug fixes to TPC files
150 Revision 1.18 2000/04/17 09:37:33 kowal2
151 removed obsolete AliTPCDigitsDisplay.C
153 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
155 New, experimental data structure from M. Ivanov
156 New tracking algorithm
157 Different pad geometry for different sectors
158 Digitization rewritten
160 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
161 Not used anymore - removed
163 Revision 1.17 2000/01/19 17:17:30 fca
164 Introducing a list of lists of hits -- more hits allowed for detector now
166 Revision 1.16 1999/11/05 09:29:23 fca
167 Accept only signals > 0
169 Revision 1.15 1999/10/08 06:26:53 fca
170 Removed ClustersIndex - not used anymore
172 Revision 1.14 1999/09/29 09:24:33 fca
173 Introduction of the Copyright and cvs Log
177 ///////////////////////////////////////////////////////////////////////////////
179 // Time Projection Chamber //
180 // This class contains the basic functions for the Time Projection Chamber //
181 // detector. Functions specific to one particular geometry are //
182 // contained in the derived classes //
186 <img src="picts/AliTPCClass.gif">
191 ///////////////////////////////////////////////////////////////////////////////
199 #include <TGeometry.h>
202 #include <TObjectTable.h>
203 #include "TParticle.h"
209 #include <iostream.h>
214 #include "AliTrackReference.h"
217 #include "AliTPCParamSR.h"
218 #include "AliTPCPRF2D.h"
219 #include "AliTPCRF1D.h"
220 #include "AliDigits.h"
221 #include "AliSimDigits.h"
222 #include "AliTPCTrackHits.h"
223 #include "AliTPCTrackHitsV2.h"
224 #include "AliPoints.h"
225 #include "AliArrayBranch.h"
228 #include "AliTPCDigitsArray.h"
229 #include "AliComplexCluster.h"
230 #include "AliClusters.h"
231 #include "AliTPCClustersRow.h"
232 #include "AliTPCClustersArray.h"
234 #include "AliTPCcluster.h"
235 #include "AliTPCclusterer.h"
236 #include "AliTPCtracker.h"
238 #include <TInterpreter.h>
245 //_____________________________________________________________________________
246 // helper class for fast matrix and vector manipulation - no range checking
247 // origin - Marian Ivanov
249 class AliTPCFastMatrix : public TMatrix {
251 AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb);
252 inline Float_t & UncheckedAt(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
253 inline Float_t UncheckedAtFast(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
256 AliTPCFastMatrix::AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb):
257 TMatrix(row_lwb, row_upb,col_lwb,col_upb)
260 //_____________________________________________________________________________
261 class AliTPCFastVector : public TVector {
263 AliTPCFastVector(Int_t size);
264 inline Float_t & UncheckedAt(Int_t index) const {return fElements[index];} //fast acces
267 AliTPCFastVector::AliTPCFastVector(Int_t size):TVector(size){
270 //_____________________________________________________________________________
274 // Default constructor
285 fHitType = 2; //default CONTAINERS - based on ROOT structure
292 //_____________________________________________________________________________
293 AliTPC::AliTPC(const char *name, const char *title)
294 : AliDetector(name,title)
297 // Standard constructor
301 // Initialise arrays of hits and digits
302 fHits = new TClonesArray("AliTPChit", 176);
303 gAlice->AddHitList(fHits);
308 fTrackHits = new AliTPCTrackHitsV2;
309 fTrackHits->SetHitPrecision(0.002);
310 fTrackHits->SetStepPrecision(0.003);
311 fTrackHits->SetMaxDistance(100);
313 fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
314 fTrackHitsOld->SetHitPrecision(0.002);
315 fTrackHitsOld->SetStepPrecision(0.003);
316 fTrackHitsOld->SetMaxDistance(100);
323 // Initialise counters
329 // Initialise color attributes
330 SetMarkerColor(kYellow);
333 // Set TPC parameters
337 if (!strcmp(title,"Default")) {
338 fTPCParam = new AliTPCParamSR;
340 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
346 //_____________________________________________________________________________
357 delete fTrackHits; //MI 15.09.2000
358 delete fTrackHitsOld; //MI 10.12.2001
359 if (fNoiseTable) delete [] fNoiseTable;
363 //_____________________________________________________________________________
364 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
367 // Add a hit to the list
369 // TClonesArray &lhits = *fHits;
370 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
372 TClonesArray &lhits = *fHits;
373 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
376 AddHit2(track,vol,hits);
379 void AliTPC::AddTrackReference(Int_t lab, TLorentzVector p, TLorentzVector x){
381 // add a trackrefernce to the list
382 if (!fTrackReferences) {
383 cerr<<"Container trackrefernce not active\n";
386 Int_t nref = fTrackReferences->GetEntriesFast();
387 TClonesArray &lref = *fTrackReferences;
388 AliTrackReference * ref = new(lref[nref]) AliTrackReference;
389 ref->SetMomentum(p[0],p[1],p[2]);
390 ref->SetPosition(x[0],x[1],x[2]);
394 //_____________________________________________________________________________
395 void AliTPC::BuildGeometry()
399 // Build TPC ROOT TNode geometry for the event display
404 const int kColorTPC=19;
405 char name[5], title[25];
406 const Double_t kDegrad=TMath::Pi()/180;
407 const Double_t kRaddeg=180./TMath::Pi();
410 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
411 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
413 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
414 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
416 Int_t nLo = fTPCParam->GetNInnerSector()/2;
417 Int_t nHi = fTPCParam->GetNOuterSector()/2;
419 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
420 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
421 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
422 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
425 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
426 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
432 // Get ALICE top node
435 nTop=gAlice->GetGeometry()->GetNode("alice");
439 rl = fTPCParam->GetInnerRadiusLow();
440 ru = fTPCParam->GetInnerRadiusUp();
444 sprintf(name,"LS%2.2d",i);
446 sprintf(title,"TPC low sector %3d",i);
449 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
450 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
451 tubs->SetNumberOfDivisions(1);
453 nNode = new TNode(name,title,name,0,0,0,"");
454 nNode->SetLineColor(kColorTPC);
460 rl = fTPCParam->GetOuterRadiusLow();
461 ru = fTPCParam->GetOuterRadiusUp();
464 sprintf(name,"US%2.2d",i);
466 sprintf(title,"TPC upper sector %d",i);
468 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
469 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
470 tubs->SetNumberOfDivisions(1);
472 nNode = new TNode(name,title,name,0,0,0,"");
473 nNode->SetLineColor(kColorTPC);
479 //_____________________________________________________________________________
480 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
483 // Calculate distance from TPC to mouse on the display
489 void AliTPC::Clusters2Tracks(TFile *of) {
490 //-----------------------------------------------------------------
491 // This is a track finder.
492 //-----------------------------------------------------------------
493 AliTPCtracker tracker(fTPCParam);
494 tracker.Clusters2Tracks(gFile,of);
497 //_____________________________________________________________________________
498 void AliTPC::CreateMaterials()
500 //-----------------------------------------------
501 // Create Materials for for TPC simulations
502 //-----------------------------------------------
504 //-----------------------------------------------------------------
505 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
506 //-----------------------------------------------------------------
508 Int_t iSXFLD=gAlice->Field()->Integ();
509 Float_t sXMGMX=gAlice->Field()->Max();
511 Float_t amat[5]; // atomic numbers
512 Float_t zmat[5]; // z
513 Float_t wmat[5]; // proportions
519 //***************** Gases *************************
521 //-------------------------------------------------
523 //-------------------------------------------------
534 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
544 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
547 //--------------------------------------------------------------
549 //--------------------------------------------------------------
566 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
568 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
583 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
585 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
601 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
603 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
605 //----------------------------------------------------------------
606 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
607 //----------------------------------------------------------------
613 Float_t rho,absl,X0,buf[1];
617 for(nc = 0;nc<fNoComp;nc++)
620 // retrive material constants
622 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
627 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
629 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
630 density += fMixtProp[nc]*rho; // density of the mixture
634 // mixture proportions by weight!
636 for(nc = 0;nc<fNoComp;nc++)
639 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
641 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
642 apure[nnc] : amol[nnc])/am;
646 // Drift gases 1 - nonsensitive, 2 - sensitive
648 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
649 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
658 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
661 //----------------------------------------------------------------------
663 //----------------------------------------------------------------------
685 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
707 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
725 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
743 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
745 // SiO2 - used later for the glass fiber
757 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
766 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
775 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
784 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
802 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
821 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
823 // Epoxy - C14 H20 O3
840 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
848 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
852 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,X0,absl,buf,nbuf);
856 wmat[0]=0.644; // by weight!
859 density=0.5*(1.25+2.265);
861 AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
865 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
867 wmat[0]=0.725; // by weight!
872 AliMixture(39,"G10",amat,zmat,density,2,wmat);
877 //----------------------------------------------------------
878 // tracking media for gases
879 //----------------------------------------------------------
881 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
882 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
883 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
884 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
886 //-----------------------------------------------------------
887 // tracking media for solids
888 //-----------------------------------------------------------
890 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
891 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
892 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
893 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
894 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
895 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
896 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
897 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
898 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
899 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
900 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
901 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
905 void AliTPC::GenerNoise(Int_t tablesize)
908 //Generate table with noise
915 if (fNoiseTable) delete[] fNoiseTable;
916 fNoiseTable = new Float_t[tablesize];
917 fNoiseDepth = tablesize;
918 fCurrentNoise =0; //!index of the noise in the noise table
920 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
921 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
924 Float_t AliTPC::GetNoise()
926 // get noise from table
927 // if ((fCurrentNoise%10)==0)
928 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
929 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
930 return fNoiseTable[fCurrentNoise++];
931 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
935 Bool_t AliTPC::IsSectorActive(Int_t sec)
938 // check if the sector is active
939 if (!fActiveSectors) return kTRUE;
940 else return fActiveSectors[sec];
943 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
945 // activate interesting sectors
946 if (fActiveSectors) delete [] fActiveSectors;
947 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
948 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
949 for (Int_t i=0;i<n;i++)
950 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
954 void AliTPC::SetActiveSectors(Int_t flag)
957 // activate sectors which were hitted by tracks
959 if (fHitType==0) return; // if Clones hit - not short volume ID information
960 if (fActiveSectors) delete [] fActiveSectors;
961 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
963 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
966 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
968 if (fHitType>1) branch = gAlice->TreeH()->GetBranch("TPC2");
969 else branch = gAlice->TreeH()->GetBranch("TPC");
970 Stat_t ntracks = gAlice->TreeH()->GetEntries();
971 // loop over all hits
972 for(Int_t track=0;track<ntracks;track++){
975 if (fTrackHits && fHitType&4) {
976 TBranch * br1 = gAlice->TreeH()->GetBranch("fVolumes");
977 TBranch * br2 = gAlice->TreeH()->GetBranch("fNVolumes");
978 br1->GetEvent(track);
979 br2->GetEvent(track);
980 Int_t *volumes = fTrackHits->GetVolumes();
981 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
982 fActiveSectors[volumes[j]]=kTRUE;
986 if (fTrackHitsOld && fHitType&2) {
987 TBranch * br = gAlice->TreeH()->GetBranch("fTrackHitsInfo");
989 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
990 for (UInt_t j=0;j<ar->GetSize();j++){
991 fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
1001 void AliTPC::Digits2Clusters(TFile *of, Int_t eventnumber)
1003 //-----------------------------------------------------------------
1004 // This is a simple cluster finder.
1005 //-----------------------------------------------------------------
1006 AliTPCclusterer::Digits2Clusters(fTPCParam,of,eventnumber);
1009 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
1010 extern Double_t SigmaZ2(Double_t, Double_t);
1011 //_____________________________________________________________________________
1012 void AliTPC::Hits2Clusters(TFile *of, Int_t eventn)
1014 //--------------------------------------------------------
1015 // TPC simple cluster generator from hits
1016 // obtained from the TPC Fast Simulator
1017 // The point errors are taken from the parametrization
1018 //--------------------------------------------------------
1020 //-----------------------------------------------------------------
1021 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1022 //-----------------------------------------------------------------
1023 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
1024 // Jouri.Belikov@cern.ch
1025 //----------------------------------------------------------------
1027 /////////////////////////////////////////////////////////////////////////////
1029 //---------------------------------------------------------------------
1030 // ALICE TPC Cluster Parameters
1031 //--------------------------------------------------------------------
1035 // Cluster width in rphi
1036 const Float_t kACrphi=0.18322;
1037 const Float_t kBCrphi=0.59551e-3;
1038 const Float_t kCCrphi=0.60952e-1;
1039 // Cluster width in z
1040 const Float_t kACz=0.19081;
1041 const Float_t kBCz=0.55938e-3;
1042 const Float_t kCCz=0.30428;
1044 TDirectory *savedir=gDirectory;
1046 if (!of->IsOpen()) {
1047 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
1051 //if(fDefaults == 0) SetDefaults();
1053 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
1055 TParticle *particle; // pointer to a given particle
1056 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1060 Float_t pl,pt,tanth,rpad,ratio;
1063 //---------------------------------------------------------------
1064 // Get the access to the tracks
1065 //---------------------------------------------------------------
1067 TTree *tH = gAlice->TreeH();
1068 Stat_t ntracks = tH->GetEntries();
1070 //Switch to the output file
1075 sprintf(cname,"TreeC_TPC_%d",eventn);
1077 fTPCParam->Write(fTPCParam->GetTitle());
1078 AliTPCClustersArray carray;
1079 carray.Setup(fTPCParam);
1080 carray.SetClusterType("AliTPCcluster");
1083 Int_t nclusters=0; //cluster counter
1085 //------------------------------------------------------------
1086 // Loop over all sectors (72 sectors for 20 deg
1087 // segmentation for both lower and upper sectors)
1088 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
1089 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1091 // First cluster for sector 0 starts at "0"
1092 //------------------------------------------------------------
1094 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
1096 fTPCParam->AdjustCosSin(isec,cph,sph);
1098 //------------------------------------------------------------
1100 //------------------------------------------------------------
1102 for(Int_t track=0;track<ntracks;track++){
1104 tH->GetEvent(track);
1106 // Get number of the TPC hits
1108 tpcHit = (AliTPChit*)FirstHit(-1);
1114 if (tpcHit->fQ == 0.) {
1115 tpcHit = (AliTPChit*) NextHit();
1116 continue; //information about track (I.Belikov)
1118 sector=tpcHit->fSector; // sector number
1121 tpcHit = (AliTPChit*) NextHit();
1124 ipart=tpcHit->Track();
1125 particle=gAlice->Particle(ipart);
1128 if(pt < 1.e-9) pt=1.e-9;
1130 tanth = TMath::Abs(tanth);
1131 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
1132 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
1134 // space-point resolutions
1136 sigmaRphi=SigmaY2(rpad,tanth,pt);
1137 sigmaZ =SigmaZ2(rpad,tanth );
1141 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
1142 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
1144 // temporary protection
1146 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
1147 if(sigmaZ < 0.) sigmaZ=0.4e-3;
1148 if(clRphi < 0.) clRphi=2.5e-3;
1149 if(clZ < 0.) clZ=2.5e-5;
1154 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
1155 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1157 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
1158 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
1159 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
1160 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
1161 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1162 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
1163 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
1164 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
1165 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
1166 xyz[2]=sigmaRphi; // fSigmaY2
1167 xyz[3]=sigmaZ; // fSigmaZ2
1168 xyz[4]=tpcHit->fQ; // q
1170 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
1171 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
1173 Int_t tracks[3]={tpcHit->Track(), -1, -1};
1174 AliTPCcluster cluster(tracks,xyz);
1176 clrow->InsertCluster(&cluster); nclusters++;
1178 tpcHit = (AliTPChit*)NextHit();
1181 } // end of loop over hits
1183 } // end of loop over tracks
1185 Int_t nrows=fTPCParam->GetNRow(isec);
1186 for (Int_t irow=0; irow<nrows; irow++) {
1187 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
1188 if (!clrow) continue;
1189 carray.StoreRow(isec,irow);
1190 carray.ClearRow(isec,irow);
1193 } // end of loop over sectors
1195 cerr<<"Number of made clusters : "<<nclusters<<" \n";
1196 carray.GetTree()->SetName(cname);
1197 carray.GetTree()->Write();
1199 savedir->cd(); //switch back to the input file
1201 } // end of function
1203 //_________________________________________________________________
1204 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1206 //--------------------------------------------------------
1207 //calculate exact cross point of track and given pad row
1208 //resulting values are expressed in "digit" coordinata
1209 //--------------------------------------------------------
1211 //-----------------------------------------------------------------
1212 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1213 //-----------------------------------------------------------------
1215 if (fClustersArray==0){
1219 TParticle *particle; // pointer to a given particle
1220 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1221 // Int_t sector,nhits;
1223 const Int_t kcmaxhits=30000;
1224 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
1225 AliTPCFastVector & xxx = *xxxx;
1226 Int_t maxhits = kcmaxhits;
1227 //construct array for each padrow
1228 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
1229 fClustersArray->CreateRow(isec,i);
1231 //---------------------------------------------------------------
1232 // Get the access to the tracks
1233 //---------------------------------------------------------------
1235 TTree *tH = gAlice->TreeH();
1236 Stat_t ntracks = tH->GetEntries();
1237 Int_t npart = gAlice->GetNtrack();
1240 if (fHitType>1) branch = tH->GetBranch("TPC2");
1241 else branch = tH->GetBranch("TPC");
1243 //------------------------------------------------------------
1245 //------------------------------------------------------------
1247 for(Int_t track=0;track<ntracks;track++){
1248 Bool_t isInSector=kTRUE;
1250 isInSector = TrackInVolume(isec,track);
1251 if (!isInSector) continue;
1253 branch->GetEntry(track); // get next track
1255 // Get number of the TPC hits and a pointer
1259 Int_t currentIndex=0;
1260 Int_t lastrow=-1; //last writen row
1264 tpcHit = (AliTPChit*)FirstHit(-1);
1267 Int_t sector=tpcHit->fSector; // sector number
1269 tpcHit = (AliTPChit*) NextHit();
1273 ipart=tpcHit->Track();
1274 if (ipart<npart) particle=gAlice->Particle(ipart);
1278 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1279 Int_t index[3]={1,isec,0};
1280 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1281 if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
1282 if (lastrow<0) lastrow=currentrow;
1283 if (currentrow==lastrow){
1284 if ( currentIndex>=maxhits){
1286 xxx.ResizeTo(4*maxhits);
1288 xxx(currentIndex*4)=x[0];
1289 xxx(currentIndex*4+1)=x[1];
1290 xxx(currentIndex*4+2)=x[2];
1291 xxx(currentIndex*4+3)=tpcHit->fQ;
1295 if (currentIndex>2){
1307 for (Int_t index=0;index<currentIndex;index++){
1309 x=x2=x3=x4=xxx(index*4);
1317 sumy+=xxx(index*4+1);
1318 sumxy+=xxx(index*4+1)*x;
1319 sumx2y+=xxx(index*4+1)*x2;
1320 sumz+=xxx(index*4+2);
1321 sumxz+=xxx(index*4+2)*x;
1322 sumx2z+=xxx(index*4+2)*x2;
1323 sumq+=xxx(index*4+3);
1325 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1326 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1327 sumx2*(sumx*sumx3-sumx2*sumx2);
1329 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1330 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1331 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1332 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1334 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1335 sumx2*(sumx*sumx2y-sumx2*sumxy);
1336 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1337 sumx2*(sumx*sumx2z-sumx2*sumxz);
1339 if (TMath::Abs(det)<0.00001){
1340 tpcHit = (AliTPChit*)NextHit();
1344 Float_t y=detay/det+centralPad;
1345 Float_t z=detaz/det;
1346 Float_t by=detby/det; //y angle
1347 Float_t bz=detbz/det; //z angle
1348 sumy/=Float_t(currentIndex);
1349 sumz/=Float_t(currentIndex);
1351 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1353 AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
1354 // AliComplexCluster cl;
1360 cl->fTracks[0]=ipart;
1364 } //end of calculating cluster for given row
1367 tpcHit = (AliTPChit*)NextHit();
1368 } // end of loop over hits
1369 } // end of loop over tracks
1370 //write padrows to tree
1371 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1372 fClustersArray->StoreRow(isec,ii);
1373 fClustersArray->ClearRow(isec,ii);
1382 void AliTPC::SDigits2Digits2(Int_t eventnumber)
1384 //create digits from summable digits
1385 GenerNoise(500000); //create teble with noise
1388 sprintf(sname,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1389 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1391 //conect tree with sSDigits
1393 if (gAlice->GetTreeDFile()) {
1394 t = (TTree *)gAlice->GetTreeDFile()->Get(sname);
1396 t = (TTree *)gDirectory->Get(sname);
1399 cerr<<"TPC tree with sdigits not found"<<endl;
1402 AliSimDigits digarr, *dummy=&digarr;
1403 t->GetBranch("Segment")->SetAddress(&dummy);
1404 Stat_t nentries = t->GetEntries();
1406 // set zero suppression
1408 fTPCParam->SetZeroSup(2);
1410 // get zero suppression
1412 Int_t zerosup = fTPCParam->GetZeroSup();
1414 //make tree with digits
1416 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1417 arr->SetClass("AliSimDigits");
1418 arr->Setup(fTPCParam);
1419 // Note that methods arr->MakeTree have different signatures
1420 if (gAlice->GetTreeDFile()) {
1421 arr->MakeTree(gAlice->GetTreeDFile());
1423 arr->MakeTree(fDigitsFile);
1426 AliTPCParam * par =fTPCParam;
1428 //Loop over segments of the TPC
1430 for (Int_t n=0; n<nentries; n++) {
1433 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1434 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
1437 if (!IsSectorActive(sec)) continue;
1438 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1439 Int_t nrows = digrow->GetNRows();
1440 Int_t ncols = digrow->GetNCols();
1442 digrow->ExpandBuffer();
1443 digarr.ExpandBuffer();
1444 digrow->ExpandTrackBuffer();
1445 digarr.ExpandTrackBuffer();
1448 Short_t * pamp0 = digarr.GetDigits();
1449 Int_t * ptracks0 = digarr.GetTracks();
1450 Short_t * pamp1 = digrow->GetDigits();
1451 Int_t * ptracks1 = digrow->GetTracks();
1452 Int_t nelems =nrows*ncols;
1453 Int_t saturation = fTPCParam->GetADCSat();
1454 //use internal structure of the AliDigits - for speed reason
1455 //if you cahnge implementation
1456 //of the Alidigits - it must be rewriten -
1457 for (Int_t i= 0; i<nelems; i++){
1458 // Float_t q = *pamp0;
1459 //q/=16.; //conversion faktor
1460 //Float_t noise= GetNoise();
1462 //q= TMath::Nint(q);
1463 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1465 if (q>saturation) q=saturation;
1467 //if (ptracks0[0]==0)
1470 ptracks1[0]=ptracks0[0];
1471 ptracks1[nelems]=ptracks0[nelems];
1472 ptracks1[2*nelems]=ptracks0[2*nelems];
1480 arr->StoreRow(sec,row);
1481 arr->ClearRow(sec,row);
1482 // cerr<<sec<<"\t"<<row<<"\n";
1489 arr->GetTree()->SetName(dname);
1490 arr->GetTree()->AutoSave();
1493 //_________________________________________
1494 void AliTPC::Merge(TTree * intree, Int_t *mask, Int_t nin, Int_t outid)
1497 //intree - pointer to array of trees with s digits
1498 //mask - mask for each
1499 //nin - number of inputs
1500 //outid - event number of the output event
1502 //create digits from summable digits
1503 //conect tree with sSDigits
1506 AliSimDigits ** digarr = new AliSimDigits*[nin];
1507 for (Int_t i1=0;i1<nin; i1++){
1509 intree[i1].GetBranch("Segment")->SetAddress(&digarr[i1]);
1511 Stat_t nentries = intree[0].GetEntries();
1513 //make tree with digits
1515 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),outid);
1516 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1517 arr->SetClass("AliSimDigits");
1518 arr->Setup(fTPCParam);
1519 arr->MakeTree(fDigitsFile);
1521 // set zero suppression
1523 fTPCParam->SetZeroSup(2);
1525 // get zero suppression
1527 Int_t zerosup = fTPCParam->GetZeroSup();
1530 AliTPCParam * par =fTPCParam;
1532 //Loop over segments of the TPC
1533 for (Int_t n=0; n<nentries; n++) {
1535 for (Int_t i=0;i<nin; i++){
1536 //connect proper digits
1537 intree[i].GetEvent(n);
1538 digarr[i]->ExpandBuffer();
1539 digarr[i]->ExpandTrackBuffer();
1542 if (!par->AdjustSectorRow(digarr[0]->GetID(),sec,row)) {
1543 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr[0]->GetID()<<endl;
1547 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1548 Int_t nrows = digrow->GetNRows();
1549 Int_t ncols = digrow->GetNCols();
1551 digrow->ExpandBuffer();
1552 digrow->ExpandTrackBuffer();
1554 for (Int_t rows=0;rows<nrows; rows++){
1555 for (Int_t col=0;col<ncols; col++){
1557 Int_t label[1000]; // i hope no more than 300 events merged
1559 // looop over digits
1560 for (Int_t i=0;i<nin; i++){
1561 q += digarr[i]->GetDigitFast(rows,col);
1562 for (Int_t tr=0;tr<3;tr++) {
1563 Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
1565 label[labptr]=lab+mask[i]-2;
1571 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac()*16);
1573 q/=16; //conversion faktor
1578 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat();
1579 digrow->SetDigitFast((Short_t)q,rows,col);
1580 for (Int_t tr=0;tr<3;tr++){
1582 ((AliSimDigits*)digrow)->SetTrackIDFast(label[tr],rows,col,tr);
1584 ((AliSimDigits*)digrow)->SetTrackIDFast(-1,rows,col,tr);
1589 arr->StoreRow(sec,row);
1591 arr->ClearRow(sec,row);
1596 arr->GetTree()->SetName(dname);
1597 arr->GetTree()->Write();
1603 /*_________________________________________
1604 void AliTPC::SDigits2Digits(Int_t eventnumber)
1608 cerr<<"Digitizing TPC...\n";
1610 Hits2Digits(eventnumber);
1615 // char treeName[100];
1617 // sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1619 // GetDigitsArray()->GetTree()->Write(treeName);
1622 //__________________________________________________________________
1623 void AliTPC::SetDefaults(){
1626 cerr<<"Setting default parameters...\n";
1628 // Set response functions
1630 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1632 printf("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...\n");
1634 param = new AliTPCParamSR();
1637 param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1640 printf("No TPC parameters found\n");
1645 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1646 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1647 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1648 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1649 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1650 rf->SetOffset(3*param->GetZSigma());
1653 TDirectory *savedir=gDirectory;
1654 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1656 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1659 prfinner->Read("prf_07504_Gati_056068_d02");
1660 prfouter1->Read("prf_10006_Gati_047051_d03");
1661 prfouter2->Read("prf_15006_Gati_047051_d03");
1665 param->SetInnerPRF(prfinner);
1666 param->SetOuter1PRF(prfouter1);
1667 param->SetOuter2PRF(prfouter2);
1668 param->SetTimeRF(rf);
1678 //__________________________________________________________________
1679 void AliTPC::Hits2Digits(Int_t eventnumber)
1681 //----------------------------------------------------
1682 // Loop over all sectors for a single event
1683 //----------------------------------------------------
1686 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1687 GenerNoise(500000); //create teble with noise
1689 //setup TPCDigitsArray
1691 if(GetDigitsArray()) delete GetDigitsArray();
1693 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1694 arr->SetClass("AliSimDigits");
1695 arr->Setup(fTPCParam);
1696 // Note that methods arr->MakeTree have different signatures
1697 if (gAlice->GetTreeDFile()) {
1698 arr->MakeTree(gAlice->GetTreeDFile());
1700 arr->MakeTree(fDigitsFile);
1702 SetDigitsArray(arr);
1704 fDigitsSwitch=0; // standard digits
1706 cerr<<"Digitizing TPC -- normal digits...\n";
1708 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1714 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1716 GetDigitsArray()->GetTree()->SetName(treeName);
1717 GetDigitsArray()->GetTree()->AutoSave();
1724 //__________________________________________________________________
1725 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1728 //-----------------------------------------------------------
1729 // summable digits - 16 bit "ADC", no noise, no saturation
1730 //-----------------------------------------------------------
1732 //----------------------------------------------------
1733 // Loop over all sectors for a single event
1734 //----------------------------------------------------
1737 if(fDefaults == 0) SetDefaults();
1738 GenerNoise(500000); //create table with noise
1739 //setup TPCDigitsArray
1741 if(GetDigitsArray()) delete GetDigitsArray();
1743 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1744 arr->SetClass("AliSimDigits");
1745 arr->Setup(fTPCParam);
1746 // Note that methods arr->MakeTree have different signatures
1747 if (gAlice->GetTreeSFile()) {
1748 arr->MakeTree(gAlice->GetTreeSFile());
1750 arr->MakeTree(fDigitsFile);
1752 SetDigitsArray(arr);
1754 cerr<<"Digitizing TPC -- summable digits...\n";
1756 fDigitsSwitch=1; // summable digits
1758 // set zero suppression to "0"
1760 fTPCParam->SetZeroSup(0);
1762 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1769 sprintf(treeName,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1771 GetDigitsArray()->GetTree()->SetName(treeName);
1772 GetDigitsArray()->GetTree()->AutoSave();
1778 //__________________________________________________________________
1779 void AliTPC::Hits2SDigits()
1782 //-----------------------------------------------------------
1783 // summable digits - 16 bit "ADC", no noise, no saturation
1784 //-----------------------------------------------------------
1786 //----------------------------------------------------
1787 // Loop over all sectors for a single event
1788 //----------------------------------------------------
1789 //MI change - for pp run
1790 Int_t eventnumber = gAlice->GetEvNumber();
1792 if(fDefaults == 0) SetDefaults();
1793 GenerNoise(500000); //create table with noise
1795 //setup TPCDigitsArray
1797 if(GetDigitsArray()) delete GetDigitsArray();
1799 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1800 arr->SetClass("AliSimDigits");
1801 arr->Setup(fTPCParam);
1802 // Note that methods arr->MakeTree have different signatures
1803 if (gAlice->GetTreeSFile()) {
1804 arr->MakeTree(gAlice->GetTreeSFile());
1806 arr->MakeTree(fDigitsFile);
1808 SetDigitsArray(arr);
1810 cerr<<"Digitizing TPC -- summable digits...\n";
1812 // fDigitsSwitch=1; // summable digits -for the moment direct
1814 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1820 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1822 GetDigitsArray()->GetTree()->SetName(treeName);
1823 GetDigitsArray()->GetTree()->AutoSave();
1828 //_____________________________________________________________________________
1829 void AliTPC::Hits2DigitsSector(Int_t isec)
1831 //-------------------------------------------------------------------
1832 // TPC conversion from hits to digits.
1833 //-------------------------------------------------------------------
1835 //-----------------------------------------------------------------
1836 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1837 //-----------------------------------------------------------------
1839 //-------------------------------------------------------
1840 // Get the access to the track hits
1841 //-------------------------------------------------------
1843 // check if the parameters are set - important if one calls this method
1844 // directly, not from the Hits2Digits
1846 if(fDefaults == 0) SetDefaults();
1848 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1849 Stat_t ntracks = tH->GetEntries();
1853 //-------------------------------------------
1854 // Only if there are any tracks...
1855 //-------------------------------------------
1859 //printf("*** Processing sector number %d ***\n",isec);
1861 Int_t nrows =fTPCParam->GetNRow(isec);
1863 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1865 MakeSector(isec,nrows,tH,ntracks,row);
1867 //--------------------------------------------------------
1868 // Digitize this sector, row by row
1869 // row[i] is the pointer to the TObjArray of AliTPCFastVectors,
1870 // each one containing electrons accepted on this
1871 // row, assigned into tracks
1872 //--------------------------------------------------------
1876 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1878 for (i=0;i<nrows;i++){
1880 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1882 DigitizeRow(i,isec,row);
1884 fDigitsArray->StoreRow(isec,i);
1886 Int_t ndig = dig->GetDigitSize();
1887 if (gDebug > 10) printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1889 fDigitsArray->ClearRow(isec,i);
1892 } // end of the sector digitization
1894 for(i=0;i<nrows;i++){
1899 delete [] row; // delete the array of pointers to TObjArray-s
1903 } // end of Hits2DigitsSector
1906 //_____________________________________________________________________________
1907 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1909 //-----------------------------------------------------------
1910 // Single row digitization, coupling from the neighbouring
1911 // rows taken into account
1912 //-----------------------------------------------------------
1914 //-----------------------------------------------------------------
1915 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1916 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1917 //-----------------------------------------------------------------
1920 Float_t zerosup = fTPCParam->GetZeroSup();
1921 // Int_t nrows =fTPCParam->GetNRow(isec);
1922 fCurrentIndex[1]= isec;
1925 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1926 Int_t nofTbins = fTPCParam->GetMaxTBin();
1927 Int_t indexRange[4];
1929 // Integrated signal for this row
1930 // and a single track signal
1933 AliTPCFastMatrix *m1 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // integrated
1934 AliTPCFastMatrix *m2 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // single
1936 AliTPCFastMatrix &total = *m1;
1938 // Array of pointers to the label-signal list
1940 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1941 Float_t **pList = new Float_t* [nofDigits];
1945 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1949 //Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1950 //Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1953 for (Int_t row= row1;row<=row2;row++){
1954 Int_t nTracks= rows[row]->GetEntries();
1955 for (i1=0;i1<nTracks;i1++){
1956 fCurrentIndex[2]= row;
1957 fCurrentIndex[3]=irow+1;
1959 m2->Zero(); // clear single track signal matrix
1960 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1961 GetList(trackLabel,nofPads,m2,indexRange,pList);
1963 else GetSignal(rows[row],i1,0,m1,indexRange);
1969 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1971 Float_t fzerosup = zerosup+0.5;
1972 for(Int_t it=0;it<nofTbins;it++){
1973 Float_t *pq = &(total.UncheckedAt(0,it));
1974 for(Int_t ip=0;ip<nofPads;ip++){
1978 if(fDigitsSwitch == 0){
1980 if(q <=fzerosup) continue; // do not fill zeros
1982 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1987 if(q <= 0.) continue; // do not fill zeros
1988 if(q>2000.) q=2000.;
1994 // "real" signal or electronic noise (list = -1)?
1997 for(Int_t j1=0;j1<3;j1++){
1998 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
2003 <A NAME="AliDigits"></A>
2004 using of AliDigits object
2007 dig->SetDigitFast((Short_t)q,it,ip);
2008 if (fDigitsArray->IsSimulated())
2010 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
2011 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
2012 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
2016 } // end of loop over time buckets
2017 } // end of lop over pads
2020 // This row has been digitized, delete nonused stuff
2023 for(lp=0;lp<nofDigits;lp++){
2024 if(pList[lp]) delete [] pList[lp];
2033 } // end of DigitizeRow
2035 //_____________________________________________________________________________
2037 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
2038 AliTPCFastMatrix *m1, AliTPCFastMatrix *m2,Int_t *indexRange)
2041 //---------------------------------------------------------------
2042 // Calculates 2-D signal (pad,time) for a single track,
2043 // returns a pointer to the signal matrix and the track label
2044 // No digitization is performed at this level!!!
2045 //---------------------------------------------------------------
2047 //-----------------------------------------------------------------
2048 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2049 // Modified: Marian Ivanov
2050 //-----------------------------------------------------------------
2052 AliTPCFastVector *tv;
2054 tv = (AliTPCFastVector*)p1->At(ntr); // pointer to a track
2055 AliTPCFastVector &v = *tv;
2057 Float_t label = v(0);
2058 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1)-1)/2;
2060 Int_t nElectrons = (tv->GetNrows()-1)/4;
2061 indexRange[0]=9999; // min pad
2062 indexRange[1]=-1; // max pad
2063 indexRange[2]=9999; //min time
2064 indexRange[3]=-1; // max time
2066 AliTPCFastMatrix &signal = *m1;
2067 AliTPCFastMatrix &total = *m2;
2069 // Loop over all electrons
2071 for(Int_t nel=0; nel<nElectrons; nel++){
2073 Float_t aval = v(idx+4);
2074 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
2075 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
2076 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
2078 Int_t *index = fTPCParam->GetResBin(0);
2079 Float_t *weight = & (fTPCParam->GetResWeight(0));
2081 if (n>0) for (Int_t i =0; i<n; i++){
2082 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
2085 Int_t time=index[2];
2086 Float_t qweight = *(weight)*eltoadcfac;
2088 if (m1!=0) signal.UncheckedAt(pad,time)+=qweight;
2089 total.UncheckedAt(pad,time)+=qweight;
2090 if (indexRange[0]>pad) indexRange[0]=pad;
2091 if (indexRange[1]<pad) indexRange[1]=pad;
2092 if (indexRange[2]>time) indexRange[2]=time;
2093 if (indexRange[3]<time) indexRange[3]=time;
2100 } // end of loop over electrons
2102 return label; // returns track label when finished
2105 //_____________________________________________________________________________
2106 void AliTPC::GetList(Float_t label,Int_t np,AliTPCFastMatrix *m,
2107 Int_t *indexRange, Float_t **pList)
2109 //----------------------------------------------------------------------
2110 // Updates the list of tracks contributing to digits for a given row
2111 //----------------------------------------------------------------------
2113 //-----------------------------------------------------------------
2114 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2115 //-----------------------------------------------------------------
2117 AliTPCFastMatrix &signal = *m;
2119 // lop over nonzero digits
2121 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
2122 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
2125 // accept only the contribution larger than 500 electrons (1/2 s_noise)
2127 if(signal(ip,it)<0.5) continue;
2130 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
2132 if(!pList[globalIndex]){
2135 // Create new list (6 elements - 3 signals and 3 labels),
2138 pList[globalIndex] = new Float_t [6];
2142 *pList[globalIndex] = -1.;
2143 *(pList[globalIndex]+1) = -1.;
2144 *(pList[globalIndex]+2) = -1.;
2145 *(pList[globalIndex]+3) = -1.;
2146 *(pList[globalIndex]+4) = -1.;
2147 *(pList[globalIndex]+5) = -1.;
2150 *pList[globalIndex] = label;
2151 *(pList[globalIndex]+3) = signal(ip,it);
2155 // check the signal magnitude
2157 Float_t highest = *(pList[globalIndex]+3);
2158 Float_t middle = *(pList[globalIndex]+4);
2159 Float_t lowest = *(pList[globalIndex]+5);
2162 // compare the new signal with already existing list
2165 if(signal(ip,it)<lowest) continue; // neglect this track
2169 if (signal(ip,it)>highest){
2170 *(pList[globalIndex]+5) = middle;
2171 *(pList[globalIndex]+4) = highest;
2172 *(pList[globalIndex]+3) = signal(ip,it);
2174 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2175 *(pList[globalIndex]+1) = *pList[globalIndex];
2176 *pList[globalIndex] = label;
2178 else if (signal(ip,it)>middle){
2179 *(pList[globalIndex]+5) = middle;
2180 *(pList[globalIndex]+4) = signal(ip,it);
2182 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2183 *(pList[globalIndex]+1) = label;
2186 *(pList[globalIndex]+5) = signal(ip,it);
2187 *(pList[globalIndex]+2) = label;
2191 } // end of loop over pads
2192 } // end of loop over time bins
2197 //___________________________________________________________________
2198 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
2199 Stat_t ntracks,TObjArray **row)
2202 //-----------------------------------------------------------------
2203 // Prepares the sector digitization, creates the vectors of
2204 // tracks for each row of this sector. The track vector
2205 // contains the track label and the position of electrons.
2206 //-----------------------------------------------------------------
2208 //-----------------------------------------------------------------
2209 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2210 //-----------------------------------------------------------------
2212 Float_t gasgain = fTPCParam->GetGasGain();
2216 AliTPChit *tpcHit; // pointer to a sigle TPC hit
2219 if (fHitType>1) branch = TH->GetBranch("TPC2");
2220 else branch = TH->GetBranch("TPC");
2223 //----------------------------------------------
2224 // Create TObjArray-s, one for each row,
2225 // each TObjArray will store the AliTPCFastVectors
2226 // of electrons, one AliTPCFastVectors per each track.
2227 //----------------------------------------------
2229 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
2230 AliTPCFastVector **tracks = new AliTPCFastVector* [nrows+2]; //pointers to the track vectors
2232 for(i=0; i<nrows+2; i++){
2233 row[i] = new TObjArray;
2240 //--------------------------------------------------------------------
2241 // Loop over tracks, the "track" contains the full history
2242 //--------------------------------------------------------------------
2244 Int_t previousTrack,currentTrack;
2245 previousTrack = -1; // nothing to store so far!
2247 for(Int_t track=0;track<ntracks;track++){
2248 Bool_t isInSector=kTRUE;
2250 isInSector = TrackInVolume(isec,track);
2251 if (!isInSector) continue;
2253 branch->GetEntry(track); // get next track
2257 tpcHit = (AliTPChit*)FirstHit(-1);
2259 //--------------------------------------------------------------
2261 //--------------------------------------------------------------
2266 Int_t sector=tpcHit->fSector; // sector number
2268 tpcHit = (AliTPChit*) NextHit();
2272 currentTrack = tpcHit->Track(); // track number
2275 if(currentTrack != previousTrack){
2277 // store already filled fTrack
2279 for(i=0;i<nrows+2;i++){
2280 if(previousTrack != -1){
2281 if(nofElectrons[i]>0){
2282 AliTPCFastVector &v = *tracks[i];
2283 v(0) = previousTrack;
2284 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2285 row[i]->Add(tracks[i]);
2288 delete tracks[i]; // delete empty AliTPCFastVector
2294 tracks[i] = new AliTPCFastVector(481); // AliTPCFastVectors for the next fTrack
2296 } // end of loop over rows
2298 previousTrack=currentTrack; // update track label
2301 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2303 //---------------------------------------------------
2304 // Calculate the electron attachment probability
2305 //---------------------------------------------------
2308 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
2309 /fTPCParam->GetDriftV();
2311 Float_t attProb = fTPCParam->GetAttCoef()*
2312 fTPCParam->GetOxyCont()*time; // fraction!
2314 //-----------------------------------------------
2315 // Loop over electrons
2316 //-----------------------------------------------
2319 for(Int_t nel=0;nel<qI;nel++){
2320 // skip if electron lost due to the attachment
2321 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2326 // protection for the nonphysical avalanche size (10**6 maximum)
2328 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2329 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2332 TransportElectron(xyz,index);
2334 fTPCParam->GetPadRow(xyz,index);
2335 // row 0 - cross talk from the innermost row
2336 // row fNRow+1 cross talk from the outermost row
2337 rowNumber = index[2]+1;
2338 //transform position to local digit coordinates
2339 //relative to nearest pad row
2340 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2342 if (isec <fTPCParam->GetNInnerSector()) {
2343 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2344 y1 = fTPCParam->GetYInner(rowNumber);
2347 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2348 y1 = fTPCParam->GetYOuter(rowNumber);
2350 // gain inefficiency at the wires edges - linear
2353 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));
2355 nofElectrons[rowNumber]++;
2356 //----------------------------------
2357 // Expand vector if necessary
2358 //----------------------------------
2359 if(nofElectrons[rowNumber]>120){
2360 Int_t range = tracks[rowNumber]->GetNrows();
2361 if((nofElectrons[rowNumber])>(range-1)/4){
2363 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
2367 AliTPCFastVector &v = *tracks[rowNumber];
2368 Int_t idx = 4*nofElectrons[rowNumber]-3;
2369 Real_t * position = &(((AliTPCFastVector&)v).UncheckedAt(idx)); //make code faster
2370 memcpy(position,xyz,4*sizeof(Float_t));
2372 } // end of loop over electrons
2374 tpcHit = (AliTPChit*)NextHit();
2376 } // end of loop over hits
2377 } // end of loop over tracks
2380 // store remaining track (the last one) if not empty
2383 for(i=0;i<nrows+2;i++){
2384 if(nofElectrons[i]>0){
2385 AliTPCFastVector &v = *tracks[i];
2386 v(0) = previousTrack;
2387 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2388 row[i]->Add(tracks[i]);
2397 delete [] nofElectrons;
2400 } // end of MakeSector
2403 //_____________________________________________________________________________
2407 // Initialise TPC detector after definition of geometry
2412 printf("\n%s: ",ClassName());
2413 for(i=0;i<35;i++) printf("*");
2414 printf(" TPC_INIT ");
2415 for(i=0;i<35;i++) printf("*");
2416 printf("\n%s: ",ClassName());
2418 for(i=0;i<80;i++) printf("*");
2423 //_____________________________________________________________________________
2424 void AliTPC::MakeBranch(Option_t* option, const char *file)
2427 // Create Tree branches for the TPC.
2429 Int_t buffersize = 4000;
2430 char branchname[10];
2431 sprintf(branchname,"%s",GetName());
2433 AliDetector::MakeBranch(option,file);
2435 const char *d = strstr(option,"D");
2437 if (fDigits && gAlice->TreeD() && d) {
2438 MakeBranchInTree(gAlice->TreeD(),
2439 branchname, &fDigits, buffersize, file);
2442 if (fHitType>1) MakeBranch2(option,file); // MI change 14.09.2000
2445 //_____________________________________________________________________________
2446 void AliTPC::ResetDigits()
2449 // Reset number of digits and the digits array for this detector
2452 if (fDigits) fDigits->Clear();
2455 //_____________________________________________________________________________
2456 void AliTPC::SetSecAL(Int_t sec)
2458 //---------------------------------------------------
2459 // Activate/deactivate selection for lower sectors
2460 //---------------------------------------------------
2462 //-----------------------------------------------------------------
2463 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2464 //-----------------------------------------------------------------
2469 //_____________________________________________________________________________
2470 void AliTPC::SetSecAU(Int_t sec)
2472 //----------------------------------------------------
2473 // Activate/deactivate selection for upper sectors
2474 //---------------------------------------------------
2476 //-----------------------------------------------------------------
2477 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2478 //-----------------------------------------------------------------
2483 //_____________________________________________________________________________
2484 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2486 //----------------------------------------
2487 // Select active lower sectors
2488 //----------------------------------------
2490 //-----------------------------------------------------------------
2491 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2492 //-----------------------------------------------------------------
2502 //_____________________________________________________________________________
2503 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2504 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2505 Int_t s11 , Int_t s12)
2507 //--------------------------------
2508 // Select active upper sectors
2509 //--------------------------------
2511 //-----------------------------------------------------------------
2512 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2513 //-----------------------------------------------------------------
2529 //_____________________________________________________________________________
2530 void AliTPC::SetSens(Int_t sens)
2533 //-------------------------------------------------------------
2534 // Activates/deactivates the sensitive strips at the center of
2535 // the pad row -- this is for the space-point resolution calculations
2536 //-------------------------------------------------------------
2538 //-----------------------------------------------------------------
2539 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2540 //-----------------------------------------------------------------
2546 void AliTPC::SetSide(Float_t side=0.)
2548 // choice of the TPC side
2553 //____________________________________________________________________________
2554 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2555 Float_t p2,Float_t p3)
2558 // gax mixture definition
2572 //_____________________________________________________________________________
2574 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2577 // electron transport taking into account:
2579 // 2.ExB at the wires
2580 // 3. nonisochronity
2582 // xyz and index must be already transformed to system 1
2585 fTPCParam->Transform1to2(xyz,index);
2588 Float_t driftl=xyz[2];
2589 if(driftl<0.01) driftl=0.01;
2590 driftl=TMath::Sqrt(driftl);
2591 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2592 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2593 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2594 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2595 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2599 if (fTPCParam->GetMWPCReadout()==kTRUE){
2600 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2601 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2603 //add nonisochronity (not implemented yet)
2606 ClassImp(AliTPCdigit)
2608 //_____________________________________________________________________________
2609 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2613 // Creates a TPC digit object
2615 fSector = digits[0];
2616 fPadRow = digits[1];
2619 fSignal = digits[4];
2625 //_____________________________________________________________________________
2626 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2630 // Creates a TPC hit object
2641 //________________________________________________________________________
2642 // Additional code because of the AliTPCTrackHitsV2
2644 void AliTPC::MakeBranch2(Option_t *option,const char *file)
2647 // Create a new branch in the current Root Tree
2648 // The branch of fHits is automatically split
2649 // MI change 14.09.2000
2650 if (fHitType<2) return;
2651 char branchname[10];
2652 sprintf(branchname,"%s2",GetName());
2654 // Get the pointer to the header
2655 const char *cH = strstr(option,"H");
2657 if (fTrackHits && gAlice->TreeH() && cH && fHitType&4) {
2658 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2659 // gAlice->TreeH(),fBufferSize,99);
2660 //TBranch * branch =
2661 gAlice->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2664 // gAlice->TreeH()->GetListOfBranches()->Add(branch);
2666 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2667 const char folder [] = "RunMC/Event/Data";
2669 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2670 Publish(folder,&fTrackHits,branchname);
2672 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2673 TDirectory *wd = gDirectory;
2675 TIter next( b->GetListOfBranches());
2676 while ((b=(TBranch*)next())) {
2681 cout << "Diverting branch " << branchname << " to file " << file << endl;
2685 if (fTrackHitsOld && gAlice->TreeH() && cH && fHitType&2) {
2686 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2687 gAlice->TreeH(),fBufferSize,99);
2688 //TBranch * branch = gAlice->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2691 gAlice->TreeH()->GetListOfBranches()->Add(branch);
2693 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2694 const char folder [] = "RunMC/Event/Data";
2696 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2697 Publish(folder,&fTrackHitsOld,branchname);
2699 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2700 TDirectory *wd = gDirectory;
2702 TIter next( b->GetListOfBranches());
2703 while ((b=(TBranch*)next())) {
2708 cout << "Diverting branch " << branchname << " to file " << file << endl;
2713 void AliTPC::SetTreeAddress()
2715 if (fHitType<=1) AliDetector::SetTreeAddress();
2716 if (fHitType>1) SetTreeAddress2();
2719 void AliTPC::SetTreeAddress2()
2722 // Set branch address for the TrackHits Tree
2725 char branchname[20];
2726 sprintf(branchname,"%s2",GetName());
2728 // Branch address for hit tree
2729 TTree *treeH = gAlice->TreeH();
2730 if ((treeH)&&(fHitType&4)) {
2731 branch = treeH->GetBranch(branchname);
2732 if (branch) branch->SetAddress(&fTrackHits);
2734 if ((treeH)&&(fHitType&2)) {
2735 branch = treeH->GetBranch(branchname);
2736 if (branch) branch->SetAddress(&fTrackHitsOld);
2741 void AliTPC::FinishPrimary()
2743 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2744 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2748 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2751 // add hit to the list
2754 int primary = gAlice->GetPrimary(track);
2755 gAlice->Particle(primary)->SetBit(kKeepBit);
2759 gAlice->FlagTrack(track);
2761 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2762 //if (hit->fTrack!=rtrack)
2763 // cout<<"bad track number\n";
2764 if (fTrackHits && fHitType&4)
2765 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2766 hits[1],hits[2],(Int_t)hits[3]);
2767 if (fTrackHitsOld &&fHitType&2 )
2768 fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2769 hits[1],hits[2],(Int_t)hits[3]);
2773 void AliTPC::ResetHits()
2775 if (fHitType&1) AliDetector::ResetHits();
2776 if (fHitType>1) ResetHits2();
2779 void AliTPC::ResetHits2()
2783 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2784 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2788 AliHit* AliTPC::FirstHit(Int_t track)
2790 if (fHitType>1) return FirstHit2(track);
2791 return AliDetector::FirstHit(track);
2793 AliHit* AliTPC::NextHit()
2795 if (fHitType>1) return NextHit2();
2797 return AliDetector::NextHit();
2800 AliHit* AliTPC::FirstHit2(Int_t track)
2803 // Initialise the hit iterator
2804 // Return the address of the first hit for track
2805 // If track>=0 the track is read from disk
2806 // while if track<0 the first hit of the current
2807 // track is returned
2810 gAlice->ResetHits();
2811 gAlice->TreeH()->GetEvent(track);
2814 if (fTrackHits && fHitType&4) {
2815 fTrackHits->First();
2816 return fTrackHits->GetHit();
2818 if (fTrackHitsOld && fHitType&2) {
2819 fTrackHitsOld->First();
2820 return fTrackHitsOld->GetHit();
2826 AliHit* AliTPC::NextHit2()
2829 //Return the next hit for the current track
2832 if (fTrackHitsOld && fHitType&2) {
2833 fTrackHitsOld->Next();
2834 return fTrackHitsOld->GetHit();
2838 return fTrackHits->GetHit();
2844 void AliTPC::LoadPoints(Int_t)
2848 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2851 if(fHitType==1) AliDetector::LoadPoints(a);
2852 else LoadPoints2(a);
2859 void AliTPC::RemapTrackHitIDs(Int_t *map)
2861 if (!fTrackHits) return;
2863 if (fTrackHitsOld && fHitType&2){
2864 AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2865 for (UInt_t i=0;i<arr->GetSize();i++){
2866 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2867 info->fTrackID = map[info->fTrackID];
2870 if (fTrackHitsOld && fHitType&4){
2871 TClonesArray * arr = fTrackHits->GetArray();;
2872 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2873 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2874 info->fTrackID = map[info->fTrackID];
2879 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2881 //return bool information - is track in given volume
2882 //load only part of the track information
2883 //return true if current track is in volume
2886 if (fTrackHitsOld && fHitType&2) {
2887 TBranch * br = gAlice->TreeH()->GetBranch("fTrackHitsInfo");
2888 br->GetEvent(track);
2889 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2890 for (UInt_t j=0;j<ar->GetSize();j++){
2891 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2895 if (fTrackHits && fHitType&4) {
2896 TBranch * br1 = gAlice->TreeH()->GetBranch("fVolumes");
2897 TBranch * br2 = gAlice->TreeH()->GetBranch("fNVolumes");
2898 br2->GetEvent(track);
2899 br1->GetEvent(track);
2900 Int_t *volumes = fTrackHits->GetVolumes();
2901 Int_t nvolumes = fTrackHits->GetNVolumes();
2902 if (!volumes && nvolumes>0) {
2903 printf("Problematic track\t%d\t%d",track,nvolumes);
2906 for (Int_t j=0;j<nvolumes; j++)
2907 if (volumes[j]==id) return kTRUE;
2911 TBranch * br = gAlice->TreeH()->GetBranch("fSector");
2912 br->GetEvent(track);
2913 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2914 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2921 //_____________________________________________________________________________
2922 void AliTPC::LoadPoints2(Int_t)
2925 // Store x, y, z of all hits in memory
2927 if (fTrackHits == 0 && fTrackHitsOld==0) return;
2930 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2931 if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2933 if (nhits == 0) return;
2934 Int_t tracks = gAlice->GetNtrack();
2935 if (fPoints == 0) fPoints = new TObjArray(tracks);
2938 Int_t *ntrk=new Int_t[tracks];
2939 Int_t *limi=new Int_t[tracks];
2940 Float_t **coor=new Float_t*[tracks];
2941 for(Int_t i=0;i<tracks;i++) {
2947 AliPoints *points = 0;
2950 Int_t chunk=nhits/4+1;
2952 // Loop over all the hits and store their position
2954 ahit = FirstHit2(-1);
2956 trk=ahit->GetTrack();
2957 if(ntrk[trk]==limi[trk]) {
2959 // Initialise a new track
2960 fp=new Float_t[3*(limi[trk]+chunk)];
2962 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2963 delete [] coor[trk];
2970 fp[3*ntrk[trk] ] = ahit->X();
2971 fp[3*ntrk[trk]+1] = ahit->Y();
2972 fp[3*ntrk[trk]+2] = ahit->Z();
2980 for(trk=0; trk<tracks; ++trk) {
2982 points = new AliPoints();
2983 points->SetMarkerColor(GetMarkerColor());
2984 points->SetMarkerSize(GetMarkerSize());
2985 points->SetDetector(this);
2986 points->SetParticle(trk);
2987 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2988 fPoints->AddAt(points,trk);
2989 delete [] coor[trk];
2999 //_____________________________________________________________________________
3000 void AliTPC::LoadPoints3(Int_t)
3003 // Store x, y, z of all hits in memory
3004 // - only intersection point with pad row
3005 if (fTrackHits == 0) return;
3007 Int_t nhits = fTrackHits->GetEntriesFast();
3008 if (nhits == 0) return;
3009 Int_t tracks = gAlice->GetNtrack();
3010 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
3011 fPoints->Expand(2*tracks);
3014 Int_t *ntrk=new Int_t[tracks];
3015 Int_t *limi=new Int_t[tracks];
3016 Float_t **coor=new Float_t*[tracks];
3017 for(Int_t i=0;i<tracks;i++) {
3023 AliPoints *points = 0;
3026 Int_t chunk=nhits/4+1;
3028 // Loop over all the hits and store their position
3030 ahit = FirstHit2(-1);
3031 //for (Int_t hit=0;hit<nhits;hit++) {
3035 // ahit = (AliHit*)fHits->UncheckedAt(hit);
3036 trk=ahit->GetTrack();
3037 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
3038 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
3039 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
3040 if (currentrow!=lastrow){
3041 lastrow = currentrow;
3042 //later calculate intersection point
3043 if(ntrk[trk]==limi[trk]) {
3045 // Initialise a new track
3046 fp=new Float_t[3*(limi[trk]+chunk)];
3048 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
3049 delete [] coor[trk];
3056 fp[3*ntrk[trk] ] = ahit->X();
3057 fp[3*ntrk[trk]+1] = ahit->Y();
3058 fp[3*ntrk[trk]+2] = ahit->Z();
3065 for(trk=0; trk<tracks; ++trk) {
3067 points = new AliPoints();
3068 points->SetMarkerColor(GetMarkerColor()+1);
3069 points->SetMarkerStyle(5);
3070 points->SetMarkerSize(0.2);
3071 points->SetDetector(this);
3072 points->SetParticle(trk);
3073 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
3074 points->SetPolyMarker(ntrk[trk],coor[trk],30);
3075 fPoints->AddAt(points,tracks+trk);
3076 delete [] coor[trk];
3087 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
3091 //fill clones array with intersection of current point with the
3096 const Int_t kcmaxhits=30000;
3097 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
3098 AliTPCFastVector & xxx = *xxxx;
3099 Int_t maxhits = kcmaxhits;
3102 AliTPChit * tpcHit=0;
3103 tpcHit = (AliTPChit*)FirstHit2(-1);
3104 Int_t currentIndex=0;
3105 Int_t lastrow=-1; //last writen row
3108 if (tpcHit==0) continue;
3109 sector=tpcHit->fSector; // sector number
3110 ipart=tpcHit->Track();
3114 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
3115 Int_t index[3]={1,sector,0};
3116 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
3117 if (currentrow<0) continue;
3118 if (lastrow<0) lastrow=currentrow;
3119 if (currentrow==lastrow){
3120 if ( currentIndex>=maxhits){
3122 xxx.ResizeTo(4*maxhits);
3124 xxx(currentIndex*4)=x[0];
3125 xxx(currentIndex*4+1)=x[1];
3126 xxx(currentIndex*4+2)=x[2];
3127 xxx(currentIndex*4+3)=tpcHit->fQ;
3131 if (currentIndex>2){
3143 for (Int_t index=0;index<currentIndex;index++){
3145 x=x2=x3=x4=xxx(index*4);
3153 sumy+=xxx(index*4+1);
3154 sumxy+=xxx(index*4+1)*x;
3155 sumx2y+=xxx(index*4+1)*x2;
3156 sumz+=xxx(index*4+2);
3157 sumxz+=xxx(index*4+2)*x;
3158 sumx2z+=xxx(index*4+2)*x2;
3159 sumq+=xxx(index*4+3);
3161 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
3162 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
3163 sumx2*(sumx*sumx3-sumx2*sumx2);
3165 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
3166 sumx2*(sumxy*sumx3-sumx2y*sumx2);
3167 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
3168 sumx2*(sumxz*sumx3-sumx2z*sumx2);
3170 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
3171 sumx2*(sumx*sumx2y-sumx2*sumxy);
3172 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
3173 sumx2*(sumx*sumx2z-sumx2*sumxz);
3175 Float_t y=detay/det+centralPad;
3176 Float_t z=detaz/det;
3177 Float_t by=detby/det; //y angle
3178 Float_t bz=detbz/det; //z angle
3179 sumy/=Float_t(currentIndex);
3180 sumz/=Float_t(currentIndex);
3182 AliComplexCluster cl;
3188 cl.fTracks[0]=ipart;
3190 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
3191 if (row!=0) row->InsertCluster(&cl);
3194 } //end of calculating cluster for given row
3196 } // end of loop over hits
3200 //_______________________________________________________________________________
3201 void AliTPC::Digits2Reco(Int_t firstevent,Int_t lastevent)
3203 // produces rec points from digits and writes them on the root file
3204 // AliTPCclusters.root
3206 TDirectory *cwd = gDirectory;
3209 AliTPCParamSR *dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60");
3211 printf("You are running 2 pad-length geom hits with 3 pad-length geom digits\n");
3213 dig = new AliTPCParamSR();
3217 dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");
3220 printf("No TPC parameters found\n");
3225 cout<<"AliTPC::Digits2Reco: TPC parameteres have been set"<<endl;
3227 if(!gSystem->Getenv("CONFIG_FILE")){
3228 out=TFile::Open("AliTPCclusters.root","recreate");
3234 tmp1=gSystem->Getenv("CONFIG_FILE_PREFIX");
3235 tmp2=gSystem->Getenv("CONFIG_OUTDIR");
3236 sprintf(tmp3,"%s%s/AliTPCclusters.root",tmp1,tmp2);
3237 out=TFile::Open(tmp3,"recreate");
3241 cout<<"AliTPC::Digits2Reco - determination of rec points begins"<<endl;
3244 for(Int_t iev=firstevent;iev<lastevent+1;iev++){
3246 printf("Processing event %d\n",iev);
3247 Digits2Clusters(out,iev);
3249 cout<<"AliTPC::Digits2Reco - determination of rec points ended"<<endl;