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.57 2002/05/07 17:23:11 kowal2
19 Linear gain inefficiency instead of the step one at the wire edges.
21 Revision 1.56 2002/04/04 16:26:33 kowal2
22 Digits (Sdigits) go to separate files now.
24 Revision 1.55 2002/03/29 06:57:45 kowal2
25 Restored backward compatibility to use the hits from Dec. 2000 production.
27 Revision 1.54 2002/03/18 17:59:13 kowal2
28 Chnges in the pad geometry - 3 pad lengths introduced.
30 Revision 1.53 2002/02/25 11:02:56 kowal2
31 Changes towards speeding up the code. Thanks to Marian Ivanov.
33 Revision 1.52 2002/02/18 09:26:09 kowal2
34 Removed compiler warning
36 Revision 1.51 2002/01/21 17:13:21 kowal2
37 New track hits using root containers. Setting active sectors added.
39 Revision 1.50 2001/12/06 14:16:19 kowal2
42 Revision 1.49 2001/11/30 11:55:37 hristov
43 Noise table created in Hits2SDigits (M.Ivanov)
45 Revision 1.48 2001/11/24 16:10:21 kowal2
48 Revision 1.47 2001/11/19 10:25:34 kowal2
49 Nearest integer instead of integer when converting to ADC counts
51 Revision 1.46 2001/11/07 06:47:12 kowal2
54 Revision 1.45 2001/11/03 13:33:48 kowal2
55 Updated algorithms in Hits2SDigits, SDigits2Digits,
59 Revision 1.44 2001/08/30 09:28:48 hristov
60 TTree names are explicitly set via SetName(name) and then Write() is called
62 Revision 1.43 2001/07/28 12:02:54 hristov
63 Branch split level set to 99
65 Revision 1.42 2001/07/28 11:38:52 hristov
66 Loop variable declared once
68 Revision 1.41 2001/07/28 10:53:50 hristov
69 Digitisation done according to the general scheme (M.Ivanov)
71 Revision 1.40 2001/07/27 13:03:14 hristov
72 Default Branch split level set to 99
74 Revision 1.39 2001/07/26 09:09:34 kowal2
75 Hits2Reco method added
77 Revision 1.38 2001/07/20 14:32:43 kowal2
78 Processing of many events possible now
80 Revision 1.37 2001/06/12 07:17:18 kowal2
81 Hits2SDigits method implemented (summable digits)
83 Revision 1.36 2001/05/16 14:57:25 alibrary
84 New files for folders and Stack
86 Revision 1.35 2001/05/08 16:02:22 kowal2
87 Updated material specifications
89 Revision 1.34 2001/05/08 15:00:15 hristov
90 Corrections for tracking in arbitrary magnenetic field. Changes towards a concept of global Alice track. Back propagation of reconstructed tracks (Yu.Belikov)
92 Revision 1.33 2001/04/03 12:40:43 kowal2
95 Revision 1.32 2001/03/12 17:47:36 hristov
96 Changes needed on Sun with CC 5.0
98 Revision 1.31 2001/03/12 08:21:50 kowal2
99 Corrected C++ bug in the material definitions
101 Revision 1.30 2001/03/01 17:34:47 kowal2
102 Correction due to the accuracy problem
104 Revision 1.29 2001/02/28 16:34:40 kowal2
105 Protection against nonphysical values of the avalanche size,
108 Revision 1.28 2001/01/26 19:57:19 hristov
109 Major upgrade of AliRoot code
111 Revision 1.27 2001/01/13 17:29:33 kowal2
112 Sun compiler correction
114 Revision 1.26 2001/01/10 07:59:43 kowal2
115 Corrections to load points from the noncompressed hits.
117 Revision 1.25 2000/11/02 07:25:31 kowal2
118 Changes due to the new hit structure.
121 Revision 1.24 2000/10/05 16:06:09 kowal2
122 Forward declarations. Changes due to a new class AliComplexCluster.
124 Revision 1.23 2000/10/02 21:28:18 fca
125 Removal of useless dependecies via forward declarations
127 Revision 1.22 2000/07/10 20:57:39 hristov
128 Update of TPC code and macros by M.Kowalski
130 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
131 Changes to obey the coding rules
133 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
134 Splitted from AliTPCtracking
136 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
137 Changed parameter settings
139 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
141 Defaults loaded automatically (hard-wired)
142 Optional parameters can be set via macro called in the constructor
144 Revision 1.19 2000/04/18 19:00:59 fca
145 Small bug fixes to TPC files
147 Revision 1.18 2000/04/17 09:37:33 kowal2
148 removed obsolete AliTPCDigitsDisplay.C
150 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
152 New, experimental data structure from M. Ivanov
153 New tracking algorithm
154 Different pad geometry for different sectors
155 Digitization rewritten
157 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
158 Not used anymore - removed
160 Revision 1.17 2000/01/19 17:17:30 fca
161 Introducing a list of lists of hits -- more hits allowed for detector now
163 Revision 1.16 1999/11/05 09:29:23 fca
164 Accept only signals > 0
166 Revision 1.15 1999/10/08 06:26:53 fca
167 Removed ClustersIndex - not used anymore
169 Revision 1.14 1999/09/29 09:24:33 fca
170 Introduction of the Copyright and cvs Log
174 ///////////////////////////////////////////////////////////////////////////////
176 // Time Projection Chamber //
177 // This class contains the basic functions for the Time Projection Chamber //
178 // detector. Functions specific to one particular geometry are //
179 // contained in the derived classes //
183 <img src="picts/AliTPCClass.gif">
188 ///////////////////////////////////////////////////////////////////////////////
196 #include <TGeometry.h>
199 #include <TObjectTable.h>
200 #include "TParticle.h"
206 #include <iostream.h>
211 #include "AliTrackReference.h"
214 #include "AliTPCParamSR.h"
215 #include "AliTPCPRF2D.h"
216 #include "AliTPCRF1D.h"
217 #include "AliDigits.h"
218 #include "AliSimDigits.h"
219 #include "AliTPCTrackHits.h"
220 #include "AliTPCTrackHitsV2.h"
221 #include "AliPoints.h"
222 #include "AliArrayBranch.h"
225 #include "AliTPCDigitsArray.h"
226 #include "AliComplexCluster.h"
227 #include "AliClusters.h"
228 #include "AliTPCClustersRow.h"
229 #include "AliTPCClustersArray.h"
231 #include "AliTPCcluster.h"
232 #include "AliTPCclusterer.h"
233 #include "AliTPCtracker.h"
235 #include <TInterpreter.h>
242 //_____________________________________________________________________________
243 // helper class for fast matrix and vector manipulation - no range checking
244 // origin - Marian Ivanov
246 class AliTPCFastMatrix : public TMatrix {
248 AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb);
249 inline Float_t & UncheckedAt(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
250 inline Float_t UncheckedAtFast(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
253 AliTPCFastMatrix::AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb):
254 TMatrix(row_lwb, row_upb,col_lwb,col_upb)
257 //_____________________________________________________________________________
258 class AliTPCFastVector : public TVector {
260 AliTPCFastVector(Int_t size);
261 inline Float_t & UncheckedAt(Int_t index) const {return fElements[index];} //fast acces
264 AliTPCFastVector::AliTPCFastVector(Int_t size):TVector(size){
267 //_____________________________________________________________________________
271 // Default constructor
282 fHitType = 2; //default CONTAINERS - based on ROOT structure
289 //_____________________________________________________________________________
290 AliTPC::AliTPC(const char *name, const char *title)
291 : AliDetector(name,title)
294 // Standard constructor
298 // Initialise arrays of hits and digits
299 fHits = new TClonesArray("AliTPChit", 176);
300 gAlice->AddHitList(fHits);
305 fTrackHits = new AliTPCTrackHitsV2;
306 fTrackHits->SetHitPrecision(0.002);
307 fTrackHits->SetStepPrecision(0.003);
308 fTrackHits->SetMaxDistance(100);
310 fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
311 fTrackHitsOld->SetHitPrecision(0.002);
312 fTrackHitsOld->SetStepPrecision(0.003);
313 fTrackHitsOld->SetMaxDistance(100);
320 // Initialise counters
326 // Initialise color attributes
327 SetMarkerColor(kYellow);
330 // Set TPC parameters
334 if (!strcmp(title,"Default")) {
335 fTPCParam = new AliTPCParamSR;
337 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
343 //_____________________________________________________________________________
354 delete fTrackHits; //MI 15.09.2000
355 delete fTrackHitsOld; //MI 10.12.2001
356 if (fNoiseTable) delete [] fNoiseTable;
360 //_____________________________________________________________________________
361 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
364 // Add a hit to the list
366 // TClonesArray &lhits = *fHits;
367 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
369 TClonesArray &lhits = *fHits;
370 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
373 AddHit2(track,vol,hits);
376 void AliTPC::AddTrackReference(Int_t lab, TLorentzVector p, TLorentzVector x){
378 // add a trackrefernce to the list
379 if (!fTrackReferences) {
380 cerr<<"Container trackrefernce not active\n";
383 Int_t nref = fTrackReferences->GetEntriesFast();
384 TClonesArray &lref = *fTrackReferences;
385 AliTrackReference * ref = new(lref[nref]) AliTrackReference;
386 ref->SetMomentum(p[0],p[1],p[2]);
387 ref->SetPosition(x[0],x[1],x[2]);
391 //_____________________________________________________________________________
392 void AliTPC::BuildGeometry()
396 // Build TPC ROOT TNode geometry for the event display
401 const int kColorTPC=19;
402 char name[5], title[25];
403 const Double_t kDegrad=TMath::Pi()/180;
404 const Double_t kRaddeg=180./TMath::Pi();
407 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
408 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
410 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
411 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
413 Int_t nLo = fTPCParam->GetNInnerSector()/2;
414 Int_t nHi = fTPCParam->GetNOuterSector()/2;
416 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
417 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
418 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
419 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
422 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
423 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
429 // Get ALICE top node
432 nTop=gAlice->GetGeometry()->GetNode("alice");
436 rl = fTPCParam->GetInnerRadiusLow();
437 ru = fTPCParam->GetInnerRadiusUp();
441 sprintf(name,"LS%2.2d",i);
443 sprintf(title,"TPC low sector %3d",i);
446 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
447 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
448 tubs->SetNumberOfDivisions(1);
450 nNode = new TNode(name,title,name,0,0,0,"");
451 nNode->SetLineColor(kColorTPC);
457 rl = fTPCParam->GetOuterRadiusLow();
458 ru = fTPCParam->GetOuterRadiusUp();
461 sprintf(name,"US%2.2d",i);
463 sprintf(title,"TPC upper sector %d",i);
465 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
466 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
467 tubs->SetNumberOfDivisions(1);
469 nNode = new TNode(name,title,name,0,0,0,"");
470 nNode->SetLineColor(kColorTPC);
476 //_____________________________________________________________________________
477 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
480 // Calculate distance from TPC to mouse on the display
486 void AliTPC::Clusters2Tracks(TFile *of) {
487 //-----------------------------------------------------------------
488 // This is a track finder.
489 //-----------------------------------------------------------------
490 AliTPCtracker tracker(fTPCParam);
491 tracker.Clusters2Tracks(gFile,of);
494 //_____________________________________________________________________________
495 void AliTPC::CreateMaterials()
497 //-----------------------------------------------
498 // Create Materials for for TPC simulations
499 //-----------------------------------------------
501 //-----------------------------------------------------------------
502 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
503 //-----------------------------------------------------------------
505 Int_t iSXFLD=gAlice->Field()->Integ();
506 Float_t sXMGMX=gAlice->Field()->Max();
508 Float_t amat[5]; // atomic numbers
509 Float_t zmat[5]; // z
510 Float_t wmat[5]; // proportions
516 //***************** Gases *************************
518 //-------------------------------------------------
520 //-------------------------------------------------
531 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
541 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
544 //--------------------------------------------------------------
546 //--------------------------------------------------------------
563 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
565 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
580 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
582 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
598 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
600 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
602 //----------------------------------------------------------------
603 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
604 //----------------------------------------------------------------
610 Float_t rho,absl,X0,buf[1];
614 for(nc = 0;nc<fNoComp;nc++)
617 // retrive material constants
619 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
624 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
626 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
627 density += fMixtProp[nc]*rho; // density of the mixture
631 // mixture proportions by weight!
633 for(nc = 0;nc<fNoComp;nc++)
636 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
638 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
639 apure[nnc] : amol[nnc])/am;
643 // Drift gases 1 - nonsensitive, 2 - sensitive
645 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
646 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
655 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
658 //----------------------------------------------------------------------
660 //----------------------------------------------------------------------
682 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
704 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
722 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
740 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
742 // SiO2 - used later for the glass fiber
754 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
763 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
772 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
781 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
799 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
818 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
820 // Epoxy - C14 H20 O3
837 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
845 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
849 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,X0,absl,buf,nbuf);
853 wmat[0]=0.644; // by weight!
856 density=0.5*(1.25+2.265);
858 AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
862 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
864 wmat[0]=0.725; // by weight!
869 AliMixture(39,"G10",amat,zmat,density,2,wmat);
874 //----------------------------------------------------------
875 // tracking media for gases
876 //----------------------------------------------------------
878 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
879 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
880 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
881 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
883 //-----------------------------------------------------------
884 // tracking media for solids
885 //-----------------------------------------------------------
887 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
888 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
889 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
890 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
891 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
892 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
893 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
894 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
895 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
896 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
897 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
898 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
902 void AliTPC::GenerNoise(Int_t tablesize)
905 //Generate table with noise
912 if (fNoiseTable) delete[] fNoiseTable;
913 fNoiseTable = new Float_t[tablesize];
914 fNoiseDepth = tablesize;
915 fCurrentNoise =0; //!index of the noise in the noise table
917 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
918 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
921 Float_t AliTPC::GetNoise()
923 // get noise from table
924 // if ((fCurrentNoise%10)==0)
925 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
926 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
927 return fNoiseTable[fCurrentNoise++];
928 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
932 Bool_t AliTPC::IsSectorActive(Int_t sec)
935 // check if the sector is active
936 if (!fActiveSectors) return kTRUE;
937 else return fActiveSectors[sec];
940 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
942 // activate interesting sectors
943 if (fActiveSectors) delete [] fActiveSectors;
944 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
945 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
946 for (Int_t i=0;i<n;i++)
947 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
951 void AliTPC::SetActiveSectors(Int_t flag)
954 // activate sectors which were hitted by tracks
956 if (fHitType==0) return; // if Clones hit - not short volume ID information
957 if (fActiveSectors) delete [] fActiveSectors;
958 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
960 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
963 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
965 if (fHitType>1) branch = gAlice->TreeH()->GetBranch("TPC2");
966 else branch = gAlice->TreeH()->GetBranch("TPC");
967 Stat_t ntracks = gAlice->TreeH()->GetEntries();
968 // loop over all hits
969 for(Int_t track=0;track<ntracks;track++){
972 if (fTrackHits && fHitType&4) {
973 TBranch * br1 = gAlice->TreeH()->GetBranch("fVolumes");
974 TBranch * br2 = gAlice->TreeH()->GetBranch("fNVolumes");
975 br1->GetEvent(track);
976 br2->GetEvent(track);
977 Int_t *volumes = fTrackHits->GetVolumes();
978 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
979 fActiveSectors[volumes[j]]=kTRUE;
983 if (fTrackHitsOld && fHitType&2) {
984 TBranch * br = gAlice->TreeH()->GetBranch("fTrackHitsInfo");
986 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
987 for (UInt_t j=0;j<ar->GetSize();j++){
988 fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
998 void AliTPC::Digits2Clusters(TFile *of, Int_t eventnumber)
1000 //-----------------------------------------------------------------
1001 // This is a simple cluster finder.
1002 //-----------------------------------------------------------------
1003 AliTPCclusterer::Digits2Clusters(fTPCParam,of,eventnumber);
1006 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
1007 extern Double_t SigmaZ2(Double_t, Double_t);
1008 //_____________________________________________________________________________
1009 void AliTPC::Hits2Clusters(TFile *of, Int_t eventn)
1011 //--------------------------------------------------------
1012 // TPC simple cluster generator from hits
1013 // obtained from the TPC Fast Simulator
1014 // The point errors are taken from the parametrization
1015 //--------------------------------------------------------
1017 //-----------------------------------------------------------------
1018 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1019 //-----------------------------------------------------------------
1020 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
1021 // Jouri.Belikov@cern.ch
1022 //----------------------------------------------------------------
1024 /////////////////////////////////////////////////////////////////////////////
1026 //---------------------------------------------------------------------
1027 // ALICE TPC Cluster Parameters
1028 //--------------------------------------------------------------------
1032 // Cluster width in rphi
1033 const Float_t kACrphi=0.18322;
1034 const Float_t kBCrphi=0.59551e-3;
1035 const Float_t kCCrphi=0.60952e-1;
1036 // Cluster width in z
1037 const Float_t kACz=0.19081;
1038 const Float_t kBCz=0.55938e-3;
1039 const Float_t kCCz=0.30428;
1041 TDirectory *savedir=gDirectory;
1043 if (!of->IsOpen()) {
1044 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
1048 //if(fDefaults == 0) SetDefaults();
1050 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
1052 TParticle *particle; // pointer to a given particle
1053 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1057 Float_t pl,pt,tanth,rpad,ratio;
1060 //---------------------------------------------------------------
1061 // Get the access to the tracks
1062 //---------------------------------------------------------------
1064 TTree *tH = gAlice->TreeH();
1065 Stat_t ntracks = tH->GetEntries();
1067 //Switch to the output file
1072 sprintf(cname,"TreeC_TPC_%d",eventn);
1074 fTPCParam->Write(fTPCParam->GetTitle());
1075 AliTPCClustersArray carray;
1076 carray.Setup(fTPCParam);
1077 carray.SetClusterType("AliTPCcluster");
1080 Int_t nclusters=0; //cluster counter
1082 //------------------------------------------------------------
1083 // Loop over all sectors (72 sectors for 20 deg
1084 // segmentation for both lower and upper sectors)
1085 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
1086 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1088 // First cluster for sector 0 starts at "0"
1089 //------------------------------------------------------------
1091 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
1093 fTPCParam->AdjustCosSin(isec,cph,sph);
1095 //------------------------------------------------------------
1097 //------------------------------------------------------------
1099 for(Int_t track=0;track<ntracks;track++){
1101 tH->GetEvent(track);
1103 // Get number of the TPC hits
1105 tpcHit = (AliTPChit*)FirstHit(-1);
1111 if (tpcHit->fQ == 0.) {
1112 tpcHit = (AliTPChit*) NextHit();
1113 continue; //information about track (I.Belikov)
1115 sector=tpcHit->fSector; // sector number
1118 tpcHit = (AliTPChit*) NextHit();
1121 ipart=tpcHit->Track();
1122 particle=gAlice->Particle(ipart);
1125 if(pt < 1.e-9) pt=1.e-9;
1127 tanth = TMath::Abs(tanth);
1128 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
1129 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
1131 // space-point resolutions
1133 sigmaRphi=SigmaY2(rpad,tanth,pt);
1134 sigmaZ =SigmaZ2(rpad,tanth );
1138 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
1139 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
1141 // temporary protection
1143 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
1144 if(sigmaZ < 0.) sigmaZ=0.4e-3;
1145 if(clRphi < 0.) clRphi=2.5e-3;
1146 if(clZ < 0.) clZ=2.5e-5;
1151 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
1152 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1154 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
1155 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
1156 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
1157 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
1158 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1159 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
1160 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
1161 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
1162 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
1163 xyz[2]=sigmaRphi; // fSigmaY2
1164 xyz[3]=sigmaZ; // fSigmaZ2
1165 xyz[4]=tpcHit->fQ; // q
1167 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
1168 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
1170 Int_t tracks[3]={tpcHit->Track(), -1, -1};
1171 AliTPCcluster cluster(tracks,xyz);
1173 clrow->InsertCluster(&cluster); nclusters++;
1175 tpcHit = (AliTPChit*)NextHit();
1178 } // end of loop over hits
1180 } // end of loop over tracks
1182 Int_t nrows=fTPCParam->GetNRow(isec);
1183 for (Int_t irow=0; irow<nrows; irow++) {
1184 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
1185 if (!clrow) continue;
1186 carray.StoreRow(isec,irow);
1187 carray.ClearRow(isec,irow);
1190 } // end of loop over sectors
1192 cerr<<"Number of made clusters : "<<nclusters<<" \n";
1193 carray.GetTree()->SetName(cname);
1194 carray.GetTree()->Write();
1196 savedir->cd(); //switch back to the input file
1198 } // end of function
1200 //_________________________________________________________________
1201 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1203 //--------------------------------------------------------
1204 //calculate exact cross point of track and given pad row
1205 //resulting values are expressed in "digit" coordinata
1206 //--------------------------------------------------------
1208 //-----------------------------------------------------------------
1209 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1210 //-----------------------------------------------------------------
1212 if (fClustersArray==0){
1216 TParticle *particle; // pointer to a given particle
1217 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1218 // Int_t sector,nhits;
1220 const Int_t kcmaxhits=30000;
1221 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
1222 AliTPCFastVector & xxx = *xxxx;
1223 Int_t maxhits = kcmaxhits;
1224 //construct array for each padrow
1225 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
1226 fClustersArray->CreateRow(isec,i);
1228 //---------------------------------------------------------------
1229 // Get the access to the tracks
1230 //---------------------------------------------------------------
1232 TTree *tH = gAlice->TreeH();
1233 Stat_t ntracks = tH->GetEntries();
1234 Int_t npart = gAlice->GetNtrack();
1237 if (fHitType>1) branch = tH->GetBranch("TPC2");
1238 else branch = tH->GetBranch("TPC");
1240 //------------------------------------------------------------
1242 //------------------------------------------------------------
1244 for(Int_t track=0;track<ntracks;track++){
1245 Bool_t isInSector=kTRUE;
1247 isInSector = TrackInVolume(isec,track);
1248 if (!isInSector) continue;
1250 branch->GetEntry(track); // get next track
1252 // Get number of the TPC hits and a pointer
1256 Int_t currentIndex=0;
1257 Int_t lastrow=-1; //last writen row
1261 tpcHit = (AliTPChit*)FirstHit(-1);
1264 Int_t sector=tpcHit->fSector; // sector number
1266 tpcHit = (AliTPChit*) NextHit();
1270 ipart=tpcHit->Track();
1271 if (ipart<npart) particle=gAlice->Particle(ipart);
1275 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1276 Int_t index[3]={1,isec,0};
1277 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1278 if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
1279 if (lastrow<0) lastrow=currentrow;
1280 if (currentrow==lastrow){
1281 if ( currentIndex>=maxhits){
1283 xxx.ResizeTo(4*maxhits);
1285 xxx(currentIndex*4)=x[0];
1286 xxx(currentIndex*4+1)=x[1];
1287 xxx(currentIndex*4+2)=x[2];
1288 xxx(currentIndex*4+3)=tpcHit->fQ;
1292 if (currentIndex>2){
1304 for (Int_t index=0;index<currentIndex;index++){
1306 x=x2=x3=x4=xxx(index*4);
1314 sumy+=xxx(index*4+1);
1315 sumxy+=xxx(index*4+1)*x;
1316 sumx2y+=xxx(index*4+1)*x2;
1317 sumz+=xxx(index*4+2);
1318 sumxz+=xxx(index*4+2)*x;
1319 sumx2z+=xxx(index*4+2)*x2;
1320 sumq+=xxx(index*4+3);
1322 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1323 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1324 sumx2*(sumx*sumx3-sumx2*sumx2);
1326 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1327 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1328 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1329 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1331 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1332 sumx2*(sumx*sumx2y-sumx2*sumxy);
1333 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1334 sumx2*(sumx*sumx2z-sumx2*sumxz);
1336 if (TMath::Abs(det)<0.00001){
1337 tpcHit = (AliTPChit*)NextHit();
1341 Float_t y=detay/det+centralPad;
1342 Float_t z=detaz/det;
1343 Float_t by=detby/det; //y angle
1344 Float_t bz=detbz/det; //z angle
1345 sumy/=Float_t(currentIndex);
1346 sumz/=Float_t(currentIndex);
1348 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1350 AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
1351 // AliComplexCluster cl;
1357 cl->fTracks[0]=ipart;
1361 } //end of calculating cluster for given row
1364 tpcHit = (AliTPChit*)NextHit();
1365 } // end of loop over hits
1366 } // end of loop over tracks
1367 //write padrows to tree
1368 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1369 fClustersArray->StoreRow(isec,ii);
1370 fClustersArray->ClearRow(isec,ii);
1379 void AliTPC::SDigits2Digits2(Int_t eventnumber)
1381 //create digits from summable digits
1382 GenerNoise(500000); //create teble with noise
1385 sprintf(sname,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1386 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1388 //conect tree with sSDigits
1390 if (gAlice->GetTreeDFile()) {
1391 t = (TTree *)gAlice->GetTreeDFile()->Get(sname);
1393 t = (TTree *)gDirectory->Get(sname);
1396 cerr<<"TPC tree with sdigits not found"<<endl;
1399 AliSimDigits digarr, *dummy=&digarr;
1400 t->GetBranch("Segment")->SetAddress(&dummy);
1401 Stat_t nentries = t->GetEntries();
1403 // set zero suppression
1405 fTPCParam->SetZeroSup(2);
1407 // get zero suppression
1409 Int_t zerosup = fTPCParam->GetZeroSup();
1411 //make tree with digits
1413 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1414 arr->SetClass("AliSimDigits");
1415 arr->Setup(fTPCParam);
1416 // Note that methods arr->MakeTree have different signatures
1417 if (gAlice->GetTreeDFile()) {
1418 arr->MakeTree(gAlice->GetTreeDFile());
1420 arr->MakeTree(fDigitsFile);
1423 AliTPCParam * par =fTPCParam;
1425 //Loop over segments of the TPC
1427 for (Int_t n=0; n<nentries; n++) {
1430 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1431 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
1434 if (!IsSectorActive(sec)) continue;
1435 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1436 Int_t nrows = digrow->GetNRows();
1437 Int_t ncols = digrow->GetNCols();
1439 digrow->ExpandBuffer();
1440 digarr.ExpandBuffer();
1441 digrow->ExpandTrackBuffer();
1442 digarr.ExpandTrackBuffer();
1445 Short_t * pamp0 = digarr.GetDigits();
1446 Int_t * ptracks0 = digarr.GetTracks();
1447 Short_t * pamp1 = digrow->GetDigits();
1448 Int_t * ptracks1 = digrow->GetTracks();
1449 Int_t nelems =nrows*ncols;
1450 Int_t saturation = fTPCParam->GetADCSat();
1451 //use internal structure of the AliDigits - for speed reason
1452 //if you cahnge implementation
1453 //of the Alidigits - it must be rewriten -
1454 for (Int_t i= 0; i<nelems; i++){
1455 // Float_t q = *pamp0;
1456 //q/=16.; //conversion faktor
1457 //Float_t noise= GetNoise();
1459 //q= TMath::Nint(q);
1460 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1462 if (q>saturation) q=saturation;
1464 //if (ptracks0[0]==0)
1467 ptracks1[0]=ptracks0[0];
1468 ptracks1[nelems]=ptracks0[nelems];
1469 ptracks1[2*nelems]=ptracks0[2*nelems];
1477 arr->StoreRow(sec,row);
1478 arr->ClearRow(sec,row);
1479 // cerr<<sec<<"\t"<<row<<"\n";
1486 arr->GetTree()->SetName(dname);
1487 arr->GetTree()->AutoSave();
1490 //_________________________________________
1491 void AliTPC::Merge(TTree * intree, Int_t *mask, Int_t nin, Int_t outid)
1494 //intree - pointer to array of trees with s digits
1495 //mask - mask for each
1496 //nin - number of inputs
1497 //outid - event number of the output event
1499 //create digits from summable digits
1500 //conect tree with sSDigits
1503 AliSimDigits ** digarr = new AliSimDigits*[nin];
1504 for (Int_t i1=0;i1<nin; i1++){
1506 intree[i1].GetBranch("Segment")->SetAddress(&digarr[i1]);
1508 Stat_t nentries = intree[0].GetEntries();
1510 //make tree with digits
1512 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),outid);
1513 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1514 arr->SetClass("AliSimDigits");
1515 arr->Setup(fTPCParam);
1516 arr->MakeTree(fDigitsFile);
1518 // set zero suppression
1520 fTPCParam->SetZeroSup(2);
1522 // get zero suppression
1524 Int_t zerosup = fTPCParam->GetZeroSup();
1527 AliTPCParam * par =fTPCParam;
1529 //Loop over segments of the TPC
1530 for (Int_t n=0; n<nentries; n++) {
1532 for (Int_t i=0;i<nin; i++){
1533 //connect proper digits
1534 intree[i].GetEvent(n);
1535 digarr[i]->ExpandBuffer();
1536 digarr[i]->ExpandTrackBuffer();
1539 if (!par->AdjustSectorRow(digarr[0]->GetID(),sec,row)) {
1540 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr[0]->GetID()<<endl;
1544 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1545 Int_t nrows = digrow->GetNRows();
1546 Int_t ncols = digrow->GetNCols();
1548 digrow->ExpandBuffer();
1549 digrow->ExpandTrackBuffer();
1551 for (Int_t rows=0;rows<nrows; rows++){
1552 for (Int_t col=0;col<ncols; col++){
1554 Int_t label[1000]; // i hope no more than 300 events merged
1556 // looop over digits
1557 for (Int_t i=0;i<nin; i++){
1558 q += digarr[i]->GetDigitFast(rows,col);
1559 for (Int_t tr=0;tr<3;tr++) {
1560 Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
1562 label[labptr]=lab+mask[i]-2;
1568 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac()*16);
1570 q/=16; //conversion faktor
1575 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat();
1576 digrow->SetDigitFast((Short_t)q,rows,col);
1577 for (Int_t tr=0;tr<3;tr++){
1579 ((AliSimDigits*)digrow)->SetTrackIDFast(label[tr],rows,col,tr);
1581 ((AliSimDigits*)digrow)->SetTrackIDFast(-1,rows,col,tr);
1586 arr->StoreRow(sec,row);
1588 arr->ClearRow(sec,row);
1593 arr->GetTree()->SetName(dname);
1594 arr->GetTree()->Write();
1600 /*_________________________________________
1601 void AliTPC::SDigits2Digits(Int_t eventnumber)
1605 cerr<<"Digitizing TPC...\n";
1607 Hits2Digits(eventnumber);
1612 // char treeName[100];
1614 // sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1616 // GetDigitsArray()->GetTree()->Write(treeName);
1619 //__________________________________________________________________
1620 void AliTPC::SetDefaults(){
1623 cerr<<"Setting default parameters...\n";
1625 // Set response functions
1627 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1629 printf("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...\n");
1631 param = new AliTPCParamSR();
1634 param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1637 printf("No TPC parameters found\n");
1642 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1643 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1644 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1645 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1646 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1647 rf->SetOffset(3*param->GetZSigma());
1650 TDirectory *savedir=gDirectory;
1651 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1653 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1656 prfinner->Read("prf_07504_Gati_056068_d02");
1657 prfouter1->Read("prf_10006_Gati_047051_d03");
1658 prfouter2->Read("prf_15006_Gati_047051_d03");
1662 param->SetInnerPRF(prfinner);
1663 param->SetOuter1PRF(prfouter1);
1664 param->SetOuter2PRF(prfouter2);
1665 param->SetTimeRF(rf);
1675 //__________________________________________________________________
1676 void AliTPC::Hits2Digits(Int_t eventnumber)
1678 //----------------------------------------------------
1679 // Loop over all sectors for a single event
1680 //----------------------------------------------------
1683 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1684 GenerNoise(500000); //create teble with noise
1686 //setup TPCDigitsArray
1688 if(GetDigitsArray()) delete GetDigitsArray();
1690 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1691 arr->SetClass("AliSimDigits");
1692 arr->Setup(fTPCParam);
1693 // Note that methods arr->MakeTree have different signatures
1694 if (gAlice->GetTreeDFile()) {
1695 arr->MakeTree(gAlice->GetTreeDFile());
1697 arr->MakeTree(fDigitsFile);
1699 SetDigitsArray(arr);
1701 fDigitsSwitch=0; // standard digits
1703 cerr<<"Digitizing TPC -- normal digits...\n";
1705 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1711 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1713 GetDigitsArray()->GetTree()->SetName(treeName);
1714 GetDigitsArray()->GetTree()->AutoSave();
1721 //__________________________________________________________________
1722 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1725 //-----------------------------------------------------------
1726 // summable digits - 16 bit "ADC", no noise, no saturation
1727 //-----------------------------------------------------------
1729 //----------------------------------------------------
1730 // Loop over all sectors for a single event
1731 //----------------------------------------------------
1734 if(fDefaults == 0) SetDefaults();
1735 GenerNoise(500000); //create table with noise
1736 //setup TPCDigitsArray
1738 if(GetDigitsArray()) delete GetDigitsArray();
1740 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1741 arr->SetClass("AliSimDigits");
1742 arr->Setup(fTPCParam);
1743 // Note that methods arr->MakeTree have different signatures
1744 if (gAlice->GetTreeSFile()) {
1745 arr->MakeTree(gAlice->GetTreeSFile());
1747 arr->MakeTree(fDigitsFile);
1749 SetDigitsArray(arr);
1751 cerr<<"Digitizing TPC -- summable digits...\n";
1753 fDigitsSwitch=1; // summable digits
1755 // set zero suppression to "0"
1757 fTPCParam->SetZeroSup(0);
1759 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1766 sprintf(treeName,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1768 GetDigitsArray()->GetTree()->SetName(treeName);
1769 GetDigitsArray()->GetTree()->AutoSave();
1775 //__________________________________________________________________
1776 void AliTPC::Hits2SDigits()
1779 //-----------------------------------------------------------
1780 // summable digits - 16 bit "ADC", no noise, no saturation
1781 //-----------------------------------------------------------
1783 //----------------------------------------------------
1784 // Loop over all sectors for a single event
1785 //----------------------------------------------------
1786 //MI change - for pp run
1787 Int_t eventnumber = gAlice->GetEvNumber();
1789 if(fDefaults == 0) SetDefaults();
1790 GenerNoise(500000); //create table with noise
1792 //setup TPCDigitsArray
1794 if(GetDigitsArray()) delete GetDigitsArray();
1796 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1797 arr->SetClass("AliSimDigits");
1798 arr->Setup(fTPCParam);
1799 // Note that methods arr->MakeTree have different signatures
1800 if (gAlice->GetTreeSFile()) {
1801 arr->MakeTree(gAlice->GetTreeSFile());
1803 arr->MakeTree(fDigitsFile);
1805 SetDigitsArray(arr);
1807 cerr<<"Digitizing TPC -- summable digits...\n";
1809 // fDigitsSwitch=1; // summable digits -for the moment direct
1811 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1817 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1819 GetDigitsArray()->GetTree()->SetName(treeName);
1820 GetDigitsArray()->GetTree()->AutoSave();
1825 //_____________________________________________________________________________
1826 void AliTPC::Hits2DigitsSector(Int_t isec)
1828 //-------------------------------------------------------------------
1829 // TPC conversion from hits to digits.
1830 //-------------------------------------------------------------------
1832 //-----------------------------------------------------------------
1833 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1834 //-----------------------------------------------------------------
1836 //-------------------------------------------------------
1837 // Get the access to the track hits
1838 //-------------------------------------------------------
1840 // check if the parameters are set - important if one calls this method
1841 // directly, not from the Hits2Digits
1843 if(fDefaults == 0) SetDefaults();
1845 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1846 Stat_t ntracks = tH->GetEntries();
1850 //-------------------------------------------
1851 // Only if there are any tracks...
1852 //-------------------------------------------
1856 //printf("*** Processing sector number %d ***\n",isec);
1858 Int_t nrows =fTPCParam->GetNRow(isec);
1860 row= new TObjArray* [nrows];
1862 MakeSector(isec,nrows,tH,ntracks,row);
1864 //--------------------------------------------------------
1865 // Digitize this sector, row by row
1866 // row[i] is the pointer to the TObjArray of AliTPCFastVectors,
1867 // each one containing electrons accepted on this
1868 // row, assigned into tracks
1869 //--------------------------------------------------------
1873 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1875 for (i=0;i<nrows;i++){
1877 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1879 DigitizeRow(i,isec,row);
1881 fDigitsArray->StoreRow(isec,i);
1883 Int_t ndig = dig->GetDigitSize();
1884 if (gDebug > 10) printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1886 fDigitsArray->ClearRow(isec,i);
1889 } // end of the sector digitization
1891 for(i=0;i<nrows;i++){
1896 delete [] row; // delete the array of pointers to TObjArray-s
1900 } // end of Hits2DigitsSector
1903 //_____________________________________________________________________________
1904 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1906 //-----------------------------------------------------------
1907 // Single row digitization, coupling from the neighbouring
1908 // rows taken into account
1909 //-----------------------------------------------------------
1911 //-----------------------------------------------------------------
1912 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1913 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1914 //-----------------------------------------------------------------
1917 Float_t zerosup = fTPCParam->GetZeroSup();
1918 Int_t nrows =fTPCParam->GetNRow(isec);
1919 fCurrentIndex[1]= isec;
1922 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1923 Int_t nofTbins = fTPCParam->GetMaxTBin();
1924 Int_t indexRange[4];
1926 // Integrated signal for this row
1927 // and a single track signal
1930 AliTPCFastMatrix *m1 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // integrated
1931 AliTPCFastMatrix *m2 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // single
1933 AliTPCFastMatrix &total = *m1;
1935 // Array of pointers to the label-signal list
1937 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1938 Float_t **pList = new Float_t* [nofDigits];
1942 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1946 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1947 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1948 for (Int_t row= row1;row<=row2;row++){
1949 Int_t nTracks= rows[row]->GetEntries();
1950 for (i1=0;i1<nTracks;i1++){
1951 fCurrentIndex[2]= row;
1952 fCurrentIndex[3]=irow;
1954 m2->Zero(); // clear single track signal matrix
1955 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1956 GetList(trackLabel,nofPads,m2,indexRange,pList);
1958 else GetSignal(rows[row],i1,0,m1,indexRange);
1964 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1966 Float_t fzerosup = zerosup+0.5;
1967 for(Int_t it=0;it<nofTbins;it++){
1968 Float_t *pq = &(total.UncheckedAt(0,it));
1969 for(Int_t ip=0;ip<nofPads;ip++){
1973 if(fDigitsSwitch == 0){
1975 if(q <=fzerosup) continue; // do not fill zeros
1977 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1982 if(q <= 0.) continue; // do not fill zeros
1983 if(q>2000.) q=2000.;
1989 // "real" signal or electronic noise (list = -1)?
1992 for(Int_t j1=0;j1<3;j1++){
1993 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1998 <A NAME="AliDigits"></A>
1999 using of AliDigits object
2002 dig->SetDigitFast((Short_t)q,it,ip);
2003 if (fDigitsArray->IsSimulated())
2005 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
2006 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
2007 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
2011 } // end of loop over time buckets
2012 } // end of lop over pads
2015 // This row has been digitized, delete nonused stuff
2018 for(lp=0;lp<nofDigits;lp++){
2019 if(pList[lp]) delete [] pList[lp];
2028 } // end of DigitizeRow
2030 //_____________________________________________________________________________
2032 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
2033 AliTPCFastMatrix *m1, AliTPCFastMatrix *m2,Int_t *indexRange)
2036 //---------------------------------------------------------------
2037 // Calculates 2-D signal (pad,time) for a single track,
2038 // returns a pointer to the signal matrix and the track label
2039 // No digitization is performed at this level!!!
2040 //---------------------------------------------------------------
2042 //-----------------------------------------------------------------
2043 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2044 // Modified: Marian Ivanov
2045 //-----------------------------------------------------------------
2047 AliTPCFastVector *tv;
2049 tv = (AliTPCFastVector*)p1->At(ntr); // pointer to a track
2050 AliTPCFastVector &v = *tv;
2052 Float_t label = v(0);
2053 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
2055 Int_t nElectrons = (tv->GetNrows()-1)/4;
2056 indexRange[0]=9999; // min pad
2057 indexRange[1]=-1; // max pad
2058 indexRange[2]=9999; //min time
2059 indexRange[3]=-1; // max time
2061 AliTPCFastMatrix &signal = *m1;
2062 AliTPCFastMatrix &total = *m2;
2064 // Loop over all electrons
2066 for(Int_t nel=0; nel<nElectrons; nel++){
2068 Float_t aval = v(idx+4);
2069 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
2070 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
2071 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
2073 Int_t *index = fTPCParam->GetResBin(0);
2074 Float_t *weight = & (fTPCParam->GetResWeight(0));
2076 if (n>0) for (Int_t i =0; i<n; i++){
2077 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
2080 Int_t time=index[2];
2081 Float_t qweight = *(weight)*eltoadcfac;
2083 if (m1!=0) signal.UncheckedAt(pad,time)+=qweight;
2084 total.UncheckedAt(pad,time)+=qweight;
2085 if (indexRange[0]>pad) indexRange[0]=pad;
2086 if (indexRange[1]<pad) indexRange[1]=pad;
2087 if (indexRange[2]>time) indexRange[2]=time;
2088 if (indexRange[3]<time) indexRange[3]=time;
2095 } // end of loop over electrons
2097 return label; // returns track label when finished
2100 //_____________________________________________________________________________
2101 void AliTPC::GetList(Float_t label,Int_t np,AliTPCFastMatrix *m,
2102 Int_t *indexRange, Float_t **pList)
2104 //----------------------------------------------------------------------
2105 // Updates the list of tracks contributing to digits for a given row
2106 //----------------------------------------------------------------------
2108 //-----------------------------------------------------------------
2109 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2110 //-----------------------------------------------------------------
2112 AliTPCFastMatrix &signal = *m;
2114 // lop over nonzero digits
2116 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
2117 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
2120 // accept only the contribution larger than 500 electrons (1/2 s_noise)
2122 if(signal(ip,it)<0.5) continue;
2125 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
2127 if(!pList[globalIndex]){
2130 // Create new list (6 elements - 3 signals and 3 labels),
2133 pList[globalIndex] = new Float_t [6];
2137 *pList[globalIndex] = -1.;
2138 *(pList[globalIndex]+1) = -1.;
2139 *(pList[globalIndex]+2) = -1.;
2140 *(pList[globalIndex]+3) = -1.;
2141 *(pList[globalIndex]+4) = -1.;
2142 *(pList[globalIndex]+5) = -1.;
2145 *pList[globalIndex] = label;
2146 *(pList[globalIndex]+3) = signal(ip,it);
2150 // check the signal magnitude
2152 Float_t highest = *(pList[globalIndex]+3);
2153 Float_t middle = *(pList[globalIndex]+4);
2154 Float_t lowest = *(pList[globalIndex]+5);
2157 // compare the new signal with already existing list
2160 if(signal(ip,it)<lowest) continue; // neglect this track
2164 if (signal(ip,it)>highest){
2165 *(pList[globalIndex]+5) = middle;
2166 *(pList[globalIndex]+4) = highest;
2167 *(pList[globalIndex]+3) = signal(ip,it);
2169 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2170 *(pList[globalIndex]+1) = *pList[globalIndex];
2171 *pList[globalIndex] = label;
2173 else if (signal(ip,it)>middle){
2174 *(pList[globalIndex]+5) = middle;
2175 *(pList[globalIndex]+4) = signal(ip,it);
2177 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2178 *(pList[globalIndex]+1) = label;
2181 *(pList[globalIndex]+5) = signal(ip,it);
2182 *(pList[globalIndex]+2) = label;
2186 } // end of loop over pads
2187 } // end of loop over time bins
2192 //___________________________________________________________________
2193 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
2194 Stat_t ntracks,TObjArray **row)
2197 //-----------------------------------------------------------------
2198 // Prepares the sector digitization, creates the vectors of
2199 // tracks for each row of this sector. The track vector
2200 // contains the track label and the position of electrons.
2201 //-----------------------------------------------------------------
2203 //-----------------------------------------------------------------
2204 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2205 //-----------------------------------------------------------------
2207 Float_t gasgain = fTPCParam->GetGasGain();
2211 AliTPChit *tpcHit; // pointer to a sigle TPC hit
2214 if (fHitType>1) branch = TH->GetBranch("TPC2");
2215 else branch = TH->GetBranch("TPC");
2218 //----------------------------------------------
2219 // Create TObjArray-s, one for each row,
2220 // each TObjArray will store the AliTPCFastVectors
2221 // of electrons, one AliTPCFastVectors per each track.
2222 //----------------------------------------------
2224 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
2225 AliTPCFastVector **tracks = new AliTPCFastVector* [nrows]; //pointers to the track vectors
2227 for(i=0; i<nrows; i++){
2228 row[i] = new TObjArray;
2235 //--------------------------------------------------------------------
2236 // Loop over tracks, the "track" contains the full history
2237 //--------------------------------------------------------------------
2239 Int_t previousTrack,currentTrack;
2240 previousTrack = -1; // nothing to store so far!
2242 for(Int_t track=0;track<ntracks;track++){
2243 Bool_t isInSector=kTRUE;
2245 isInSector = TrackInVolume(isec,track);
2246 if (!isInSector) continue;
2248 branch->GetEntry(track); // get next track
2252 tpcHit = (AliTPChit*)FirstHit(-1);
2254 //--------------------------------------------------------------
2256 //--------------------------------------------------------------
2261 Int_t sector=tpcHit->fSector; // sector number
2263 tpcHit = (AliTPChit*) NextHit();
2267 currentTrack = tpcHit->Track(); // track number
2270 if(currentTrack != previousTrack){
2272 // store already filled fTrack
2274 for(i=0;i<nrows;i++){
2275 if(previousTrack != -1){
2276 if(nofElectrons[i]>0){
2277 AliTPCFastVector &v = *tracks[i];
2278 v(0) = previousTrack;
2279 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2280 row[i]->Add(tracks[i]);
2283 delete tracks[i]; // delete empty AliTPCFastVector
2289 tracks[i] = new AliTPCFastVector(481); // AliTPCFastVectors for the next fTrack
2291 } // end of loop over rows
2293 previousTrack=currentTrack; // update track label
2296 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2298 //---------------------------------------------------
2299 // Calculate the electron attachment probability
2300 //---------------------------------------------------
2303 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
2304 /fTPCParam->GetDriftV();
2306 Float_t attProb = fTPCParam->GetAttCoef()*
2307 fTPCParam->GetOxyCont()*time; // fraction!
2309 //-----------------------------------------------
2310 // Loop over electrons
2311 //-----------------------------------------------
2314 for(Int_t nel=0;nel<qI;nel++){
2315 // skip if electron lost due to the attachment
2316 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2321 // protection for the nonphysical avalanche size (10**6 maximum)
2323 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2324 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2327 TransportElectron(xyz,index); //MI change -august
2329 fTPCParam->GetPadRow(xyz,index); //MI change august
2330 rowNumber = index[2];
2331 //transform position to local digit coordinates
2332 //relative to nearest pad row
2333 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
2335 if (isec <fTPCParam->GetNInnerSector()) {
2336 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2337 y1 = fTPCParam->GetYInner(rowNumber);
2340 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2341 y1 = fTPCParam->GetYOuter(rowNumber);
2344 // gain inefficiency at the wires edges - linear
2348 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));
2350 nofElectrons[rowNumber]++;
2351 //----------------------------------
2352 // Expand vector if necessary
2353 //----------------------------------
2354 if(nofElectrons[rowNumber]>120){
2355 Int_t range = tracks[rowNumber]->GetNrows();
2356 if((nofElectrons[rowNumber])>(range-1)/4){
2358 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
2362 AliTPCFastVector &v = *tracks[rowNumber];
2363 Int_t idx = 4*nofElectrons[rowNumber]-3;
2364 Real_t * position = &(((AliTPCFastVector&)v).UncheckedAt(idx)); //make code faster
2365 memcpy(position,xyz,4*sizeof(Float_t));
2367 } // end of loop over electrons
2369 tpcHit = (AliTPChit*)NextHit();
2371 } // end of loop over hits
2372 } // end of loop over tracks
2375 // store remaining track (the last one) if not empty
2378 for(i=0;i<nrows;i++){
2379 if(nofElectrons[i]>0){
2380 AliTPCFastVector &v = *tracks[i];
2381 v(0) = previousTrack;
2382 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2383 row[i]->Add(tracks[i]);
2392 delete [] nofElectrons;
2395 } // end of MakeSector
2398 //_____________________________________________________________________________
2402 // Initialise TPC detector after definition of geometry
2407 printf("\n%s: ",ClassName());
2408 for(i=0;i<35;i++) printf("*");
2409 printf(" TPC_INIT ");
2410 for(i=0;i<35;i++) printf("*");
2411 printf("\n%s: ",ClassName());
2413 for(i=0;i<80;i++) printf("*");
2418 //_____________________________________________________________________________
2419 void AliTPC::MakeBranch(Option_t* option, const char *file)
2422 // Create Tree branches for the TPC.
2424 Int_t buffersize = 4000;
2425 char branchname[10];
2426 sprintf(branchname,"%s",GetName());
2428 AliDetector::MakeBranch(option,file);
2430 const char *d = strstr(option,"D");
2432 if (fDigits && gAlice->TreeD() && d) {
2433 MakeBranchInTree(gAlice->TreeD(),
2434 branchname, &fDigits, buffersize, file);
2437 if (fHitType>1) MakeBranch2(option,file); // MI change 14.09.2000
2440 //_____________________________________________________________________________
2441 void AliTPC::ResetDigits()
2444 // Reset number of digits and the digits array for this detector
2447 if (fDigits) fDigits->Clear();
2450 //_____________________________________________________________________________
2451 void AliTPC::SetSecAL(Int_t sec)
2453 //---------------------------------------------------
2454 // Activate/deactivate selection for lower sectors
2455 //---------------------------------------------------
2457 //-----------------------------------------------------------------
2458 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2459 //-----------------------------------------------------------------
2464 //_____________________________________________________________________________
2465 void AliTPC::SetSecAU(Int_t sec)
2467 //----------------------------------------------------
2468 // Activate/deactivate selection for upper sectors
2469 //---------------------------------------------------
2471 //-----------------------------------------------------------------
2472 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2473 //-----------------------------------------------------------------
2478 //_____________________________________________________________________________
2479 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2481 //----------------------------------------
2482 // Select active lower sectors
2483 //----------------------------------------
2485 //-----------------------------------------------------------------
2486 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2487 //-----------------------------------------------------------------
2497 //_____________________________________________________________________________
2498 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2499 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2500 Int_t s11 , Int_t s12)
2502 //--------------------------------
2503 // Select active upper sectors
2504 //--------------------------------
2506 //-----------------------------------------------------------------
2507 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2508 //-----------------------------------------------------------------
2524 //_____________________________________________________________________________
2525 void AliTPC::SetSens(Int_t sens)
2528 //-------------------------------------------------------------
2529 // Activates/deactivates the sensitive strips at the center of
2530 // the pad row -- this is for the space-point resolution calculations
2531 //-------------------------------------------------------------
2533 //-----------------------------------------------------------------
2534 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2535 //-----------------------------------------------------------------
2541 void AliTPC::SetSide(Float_t side=0.)
2543 // choice of the TPC side
2548 //____________________________________________________________________________
2549 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2550 Float_t p2,Float_t p3)
2553 // gax mixture definition
2567 //_____________________________________________________________________________
2569 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2572 // electron transport taking into account:
2574 // 2.ExB at the wires
2575 // 3. nonisochronity
2577 // xyz and index must be already transformed to system 1
2580 fTPCParam->Transform1to2(xyz,index);
2583 Float_t driftl=xyz[2];
2584 if(driftl<0.01) driftl=0.01;
2585 driftl=TMath::Sqrt(driftl);
2586 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2587 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2588 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2589 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2590 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2594 if (fTPCParam->GetMWPCReadout()==kTRUE){
2596 fTPCParam->Transform2to2NearestWire(xyz,index);
2597 Float_t dx=xyz[0]-x1;
2598 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2600 //add nonisochronity (not implemented yet)
2604 ClassImp(AliTPCdigit)
2606 //_____________________________________________________________________________
2607 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2611 // Creates a TPC digit object
2613 fSector = digits[0];
2614 fPadRow = digits[1];
2617 fSignal = digits[4];
2623 //_____________________________________________________________________________
2624 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2628 // Creates a TPC hit object
2639 //________________________________________________________________________
2640 // Additional code because of the AliTPCTrackHitsV2
2642 void AliTPC::MakeBranch2(Option_t *option,const char *file)
2645 // Create a new branch in the current Root Tree
2646 // The branch of fHits is automatically split
2647 // MI change 14.09.2000
2648 if (fHitType<2) return;
2649 char branchname[10];
2650 sprintf(branchname,"%s2",GetName());
2652 // Get the pointer to the header
2653 const char *cH = strstr(option,"H");
2655 if (fTrackHits && gAlice->TreeH() && cH && fHitType&4) {
2656 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2657 // gAlice->TreeH(),fBufferSize,99);
2658 //TBranch * branch =
2659 gAlice->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2662 // gAlice->TreeH()->GetListOfBranches()->Add(branch);
2664 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2665 const char folder [] = "RunMC/Event/Data";
2667 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2668 Publish(folder,&fTrackHits,branchname);
2670 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2671 TDirectory *wd = gDirectory;
2673 TIter next( b->GetListOfBranches());
2674 while ((b=(TBranch*)next())) {
2679 cout << "Diverting branch " << branchname << " to file " << file << endl;
2683 if (fTrackHitsOld && gAlice->TreeH() && cH && fHitType&2) {
2684 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2685 gAlice->TreeH(),fBufferSize,99);
2686 //TBranch * branch = gAlice->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2689 gAlice->TreeH()->GetListOfBranches()->Add(branch);
2691 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2692 const char folder [] = "RunMC/Event/Data";
2694 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2695 Publish(folder,&fTrackHitsOld,branchname);
2697 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2698 TDirectory *wd = gDirectory;
2700 TIter next( b->GetListOfBranches());
2701 while ((b=(TBranch*)next())) {
2706 cout << "Diverting branch " << branchname << " to file " << file << endl;
2711 void AliTPC::SetTreeAddress()
2713 if (fHitType<=1) AliDetector::SetTreeAddress();
2714 if (fHitType>1) SetTreeAddress2();
2717 void AliTPC::SetTreeAddress2()
2720 // Set branch address for the TrackHits Tree
2723 char branchname[20];
2724 sprintf(branchname,"%s2",GetName());
2726 // Branch address for hit tree
2727 TTree *treeH = gAlice->TreeH();
2728 if ((treeH)&&(fHitType&4)) {
2729 branch = treeH->GetBranch(branchname);
2730 if (branch) branch->SetAddress(&fTrackHits);
2732 if ((treeH)&&(fHitType&2)) {
2733 branch = treeH->GetBranch(branchname);
2734 if (branch) branch->SetAddress(&fTrackHitsOld);
2739 void AliTPC::FinishPrimary()
2741 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2742 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2746 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2749 // add hit to the list
2752 int primary = gAlice->GetPrimary(track);
2753 gAlice->Particle(primary)->SetBit(kKeepBit);
2757 gAlice->FlagTrack(track);
2759 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2760 //if (hit->fTrack!=rtrack)
2761 // cout<<"bad track number\n";
2762 if (fTrackHits && fHitType&4)
2763 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2764 hits[1],hits[2],(Int_t)hits[3]);
2765 if (fTrackHitsOld &&fHitType&2 )
2766 fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2767 hits[1],hits[2],(Int_t)hits[3]);
2771 void AliTPC::ResetHits()
2773 if (fHitType&1) AliDetector::ResetHits();
2774 if (fHitType>1) ResetHits2();
2777 void AliTPC::ResetHits2()
2781 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2782 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2786 AliHit* AliTPC::FirstHit(Int_t track)
2788 if (fHitType>1) return FirstHit2(track);
2789 return AliDetector::FirstHit(track);
2791 AliHit* AliTPC::NextHit()
2793 if (fHitType>1) return NextHit2();
2795 return AliDetector::NextHit();
2798 AliHit* AliTPC::FirstHit2(Int_t track)
2801 // Initialise the hit iterator
2802 // Return the address of the first hit for track
2803 // If track>=0 the track is read from disk
2804 // while if track<0 the first hit of the current
2805 // track is returned
2808 gAlice->ResetHits();
2809 gAlice->TreeH()->GetEvent(track);
2812 if (fTrackHits && fHitType&4) {
2813 fTrackHits->First();
2814 return fTrackHits->GetHit();
2816 if (fTrackHitsOld && fHitType&2) {
2817 fTrackHitsOld->First();
2818 return fTrackHitsOld->GetHit();
2824 AliHit* AliTPC::NextHit2()
2827 //Return the next hit for the current track
2830 if (fTrackHitsOld && fHitType&2) {
2831 fTrackHitsOld->Next();
2832 return fTrackHitsOld->GetHit();
2836 return fTrackHits->GetHit();
2842 void AliTPC::LoadPoints(Int_t)
2846 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2849 if(fHitType==1) AliDetector::LoadPoints(a);
2850 else LoadPoints2(a);
2857 void AliTPC::RemapTrackHitIDs(Int_t *map)
2859 if (!fTrackHits) return;
2861 if (fTrackHitsOld && fHitType&2){
2862 AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2863 for (UInt_t i=0;i<arr->GetSize();i++){
2864 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2865 info->fTrackID = map[info->fTrackID];
2868 if (fTrackHitsOld && fHitType&4){
2869 TClonesArray * arr = fTrackHits->GetArray();;
2870 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2871 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2872 info->fTrackID = map[info->fTrackID];
2877 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2879 //return bool information - is track in given volume
2880 //load only part of the track information
2881 //return true if current track is in volume
2884 if (fTrackHitsOld && fHitType&2) {
2885 TBranch * br = gAlice->TreeH()->GetBranch("fTrackHitsInfo");
2886 br->GetEvent(track);
2887 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2888 for (UInt_t j=0;j<ar->GetSize();j++){
2889 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2893 if (fTrackHits && fHitType&4) {
2894 TBranch * br1 = gAlice->TreeH()->GetBranch("fVolumes");
2895 TBranch * br2 = gAlice->TreeH()->GetBranch("fNVolumes");
2896 br2->GetEvent(track);
2897 br1->GetEvent(track);
2898 Int_t *volumes = fTrackHits->GetVolumes();
2899 Int_t nvolumes = fTrackHits->GetNVolumes();
2900 if (!volumes && nvolumes>0) {
2901 printf("Problematic track\t%d\t%d",track,nvolumes);
2904 for (Int_t j=0;j<nvolumes; j++)
2905 if (volumes[j]==id) return kTRUE;
2909 TBranch * br = gAlice->TreeH()->GetBranch("fSector");
2910 br->GetEvent(track);
2911 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2912 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2919 //_____________________________________________________________________________
2920 void AliTPC::LoadPoints2(Int_t)
2923 // Store x, y, z of all hits in memory
2925 if (fTrackHits == 0 && fTrackHitsOld==0) return;
2928 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2929 if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2931 if (nhits == 0) return;
2932 Int_t tracks = gAlice->GetNtrack();
2933 if (fPoints == 0) fPoints = new TObjArray(tracks);
2936 Int_t *ntrk=new Int_t[tracks];
2937 Int_t *limi=new Int_t[tracks];
2938 Float_t **coor=new Float_t*[tracks];
2939 for(Int_t i=0;i<tracks;i++) {
2945 AliPoints *points = 0;
2948 Int_t chunk=nhits/4+1;
2950 // Loop over all the hits and store their position
2952 ahit = FirstHit2(-1);
2954 trk=ahit->GetTrack();
2955 if(ntrk[trk]==limi[trk]) {
2957 // Initialise a new track
2958 fp=new Float_t[3*(limi[trk]+chunk)];
2960 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2961 delete [] coor[trk];
2968 fp[3*ntrk[trk] ] = ahit->X();
2969 fp[3*ntrk[trk]+1] = ahit->Y();
2970 fp[3*ntrk[trk]+2] = ahit->Z();
2978 for(trk=0; trk<tracks; ++trk) {
2980 points = new AliPoints();
2981 points->SetMarkerColor(GetMarkerColor());
2982 points->SetMarkerSize(GetMarkerSize());
2983 points->SetDetector(this);
2984 points->SetParticle(trk);
2985 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2986 fPoints->AddAt(points,trk);
2987 delete [] coor[trk];
2997 //_____________________________________________________________________________
2998 void AliTPC::LoadPoints3(Int_t)
3001 // Store x, y, z of all hits in memory
3002 // - only intersection point with pad row
3003 if (fTrackHits == 0) return;
3005 Int_t nhits = fTrackHits->GetEntriesFast();
3006 if (nhits == 0) return;
3007 Int_t tracks = gAlice->GetNtrack();
3008 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
3009 fPoints->Expand(2*tracks);
3012 Int_t *ntrk=new Int_t[tracks];
3013 Int_t *limi=new Int_t[tracks];
3014 Float_t **coor=new Float_t*[tracks];
3015 for(Int_t i=0;i<tracks;i++) {
3021 AliPoints *points = 0;
3024 Int_t chunk=nhits/4+1;
3026 // Loop over all the hits and store their position
3028 ahit = FirstHit2(-1);
3029 //for (Int_t hit=0;hit<nhits;hit++) {
3033 // ahit = (AliHit*)fHits->UncheckedAt(hit);
3034 trk=ahit->GetTrack();
3035 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
3036 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
3037 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
3038 if (currentrow!=lastrow){
3039 lastrow = currentrow;
3040 //later calculate intersection point
3041 if(ntrk[trk]==limi[trk]) {
3043 // Initialise a new track
3044 fp=new Float_t[3*(limi[trk]+chunk)];
3046 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
3047 delete [] coor[trk];
3054 fp[3*ntrk[trk] ] = ahit->X();
3055 fp[3*ntrk[trk]+1] = ahit->Y();
3056 fp[3*ntrk[trk]+2] = ahit->Z();
3063 for(trk=0; trk<tracks; ++trk) {
3065 points = new AliPoints();
3066 points->SetMarkerColor(GetMarkerColor()+1);
3067 points->SetMarkerStyle(5);
3068 points->SetMarkerSize(0.2);
3069 points->SetDetector(this);
3070 points->SetParticle(trk);
3071 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
3072 points->SetPolyMarker(ntrk[trk],coor[trk],30);
3073 fPoints->AddAt(points,tracks+trk);
3074 delete [] coor[trk];
3085 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
3089 //fill clones array with intersection of current point with the
3094 const Int_t kcmaxhits=30000;
3095 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
3096 AliTPCFastVector & xxx = *xxxx;
3097 Int_t maxhits = kcmaxhits;
3100 AliTPChit * tpcHit=0;
3101 tpcHit = (AliTPChit*)FirstHit2(-1);
3102 Int_t currentIndex=0;
3103 Int_t lastrow=-1; //last writen row
3106 if (tpcHit==0) continue;
3107 sector=tpcHit->fSector; // sector number
3108 ipart=tpcHit->Track();
3112 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
3113 Int_t index[3]={1,sector,0};
3114 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
3115 if (currentrow<0) continue;
3116 if (lastrow<0) lastrow=currentrow;
3117 if (currentrow==lastrow){
3118 if ( currentIndex>=maxhits){
3120 xxx.ResizeTo(4*maxhits);
3122 xxx(currentIndex*4)=x[0];
3123 xxx(currentIndex*4+1)=x[1];
3124 xxx(currentIndex*4+2)=x[2];
3125 xxx(currentIndex*4+3)=tpcHit->fQ;
3129 if (currentIndex>2){
3141 for (Int_t index=0;index<currentIndex;index++){
3143 x=x2=x3=x4=xxx(index*4);
3151 sumy+=xxx(index*4+1);
3152 sumxy+=xxx(index*4+1)*x;
3153 sumx2y+=xxx(index*4+1)*x2;
3154 sumz+=xxx(index*4+2);
3155 sumxz+=xxx(index*4+2)*x;
3156 sumx2z+=xxx(index*4+2)*x2;
3157 sumq+=xxx(index*4+3);
3159 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
3160 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
3161 sumx2*(sumx*sumx3-sumx2*sumx2);
3163 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
3164 sumx2*(sumxy*sumx3-sumx2y*sumx2);
3165 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
3166 sumx2*(sumxz*sumx3-sumx2z*sumx2);
3168 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
3169 sumx2*(sumx*sumx2y-sumx2*sumxy);
3170 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
3171 sumx2*(sumx*sumx2z-sumx2*sumxz);
3173 Float_t y=detay/det+centralPad;
3174 Float_t z=detaz/det;
3175 Float_t by=detby/det; //y angle
3176 Float_t bz=detbz/det; //z angle
3177 sumy/=Float_t(currentIndex);
3178 sumz/=Float_t(currentIndex);
3180 AliComplexCluster cl;
3186 cl.fTracks[0]=ipart;
3188 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
3189 if (row!=0) row->InsertCluster(&cl);
3192 } //end of calculating cluster for given row
3194 } // end of loop over hits
3198 //_______________________________________________________________________________
3199 void AliTPC::Digits2Reco(Int_t firstevent,Int_t lastevent)
3201 // produces rec points from digits and writes them on the root file
3202 // AliTPCclusters.root
3204 TDirectory *cwd = gDirectory;
3207 AliTPCParamSR *dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60");
3209 printf("You are running 2 pad-length geom hits with 3 pad-length geom digits\n");
3211 dig = new AliTPCParamSR();
3215 dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");
3218 printf("No TPC parameters found\n");
3223 cout<<"AliTPC::Digits2Reco: TPC parameteres have been set"<<endl;
3225 if(!gSystem->Getenv("CONFIG_FILE")){
3226 out=TFile::Open("AliTPCclusters.root","recreate");
3232 tmp1=gSystem->Getenv("CONFIG_FILE_PREFIX");
3233 tmp2=gSystem->Getenv("CONFIG_OUTDIR");
3234 sprintf(tmp3,"%s%s/AliTPCclusters.root",tmp1,tmp2);
3235 out=TFile::Open(tmp3,"recreate");
3239 cout<<"AliTPC::Digits2Reco - determination of rec points begins"<<endl;
3242 for(Int_t iev=firstevent;iev<lastevent+1;iev++){
3244 printf("Processing event %d\n",iev);
3245 Digits2Clusters(out,iev);
3247 cout<<"AliTPC::Digits2Reco - determination of rec points ended"<<endl;