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 ///////////////////////////////////////////////////////////////////////////////
20 // Time Projection Chamber //
21 // This class contains the basic functions for the Time Projection Chamber //
22 // detector. Functions specific to one particular geometry are //
23 // contained in the derived classes //
27 <img src="picts/AliTPCClass.gif">
32 ///////////////////////////////////////////////////////////////////////////////
36 #include <Riostream.h>
40 #include <TGeometry.h>
41 #include <TInterpreter.h>
46 #include <TObjectTable.h>
47 #include <TParticle.h>
53 #include <TVirtualMC.h>
56 #include <TStopwatch.h>
58 #include "AliDigits.h"
60 #include "AliPoints.h"
62 #include "AliRunLoader.h"
63 #include "AliSimDigits.h"
66 #include "AliTPCDigitsArray.h"
67 #include "AliTPCLoader.h"
68 #include "AliTPCPRF2D.h"
69 #include "AliTPCParamSR.h"
70 #include "AliTPCRF1D.h"
71 //#include "AliTPCTrackHits.h"
72 #include "AliTPCTrackHitsV2.h"
73 #include "AliTrackReference.h"
76 #include "AliTPCDigitizer.h"
77 #include "AliTPCBuffer.h"
78 #include "AliTPCDDLRawData.h"
80 #include "AliTPCcalibDB.h"
81 #include "AliTPCCalPad.h"
82 #include "AliTPCCalROC.h"
83 #include "AliTPCExB.h"
84 #include "AliRawReader.h"
85 #include "AliTPCRawStream.h"
88 //_____________________________________________________________________________
89 AliTPC::AliTPC():AliDetector(),
99 fPrimaryIonisation(0),
109 // Default constructor
113 // fTrackHitsOld = 0;
114 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
115 fHitType = 4; // ROOT containers
117 fHitType = 2; //default CONTAINERS - based on ROOT structure
121 //_____________________________________________________________________________
122 AliTPC::AliTPC(const char *name, const char *title)
123 : AliDetector(name,title),
133 fPrimaryIonisation(0),
142 // Standard constructor
146 // Initialise arrays of hits and digits
147 fHits = new TClonesArray("AliTPChit", 176);
148 gAlice->GetMCApp()->AddHitList(fHits);
150 fTrackHits = new AliTPCTrackHitsV2;
151 fTrackHits->SetHitPrecision(0.002);
152 fTrackHits->SetStepPrecision(0.003);
153 fTrackHits->SetMaxDistance(100);
155 //fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
156 //fTrackHitsOld->SetHitPrecision(0.002);
157 //fTrackHitsOld->SetStepPrecision(0.003);
158 //fTrackHitsOld->SetMaxDistance(100);
161 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
162 fHitType = 4; // ROOT containers
172 // Initialise color attributes
173 //PH SetMarkerColor(kYellow);
176 // Set TPC parameters
180 if (!strcmp(title,"Default")) {
181 fTPCParam = new AliTPCParamSR;
183 AliWarning("In Config.C you must set non-default parameters.");
188 //_____________________________________________________________________________
199 delete fTrackHits; //MI 15.09.2000
200 // delete fTrackHitsOld; //MI 10.12.2001
203 delete [] fNoiseTable;
204 delete [] fActiveSectors;
208 //_____________________________________________________________________________
209 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
212 // Add a hit to the list
215 TClonesArray &lhits = *fHits;
216 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
219 AddHit2(track,vol,hits);
222 //_____________________________________________________________________________
223 void AliTPC::BuildGeometry()
227 // Build TPC ROOT TNode geometry for the event display
232 const int kColorTPC=19;
233 char name[5], title[25];
234 const Double_t kDegrad=TMath::Pi()/180;
235 const Double_t kRaddeg=180./TMath::Pi();
238 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
239 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
241 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
242 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
244 Int_t nLo = fTPCParam->GetNInnerSector()/2;
245 Int_t nHi = fTPCParam->GetNOuterSector()/2;
247 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
248 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
249 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
250 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
253 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
254 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
260 // Get ALICE top node
263 nTop=gAlice->GetGeometry()->GetNode("alice");
267 rl = fTPCParam->GetInnerRadiusLow();
268 ru = fTPCParam->GetInnerRadiusUp();
272 sprintf(name,"LS%2.2d",i);
274 sprintf(title,"TPC low sector %3d",i);
277 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
278 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
279 tubs->SetNumberOfDivisions(1);
281 nNode = new TNode(name,title,name,0,0,0,"");
282 nNode->SetLineColor(kColorTPC);
288 rl = fTPCParam->GetOuterRadiusLow();
289 ru = fTPCParam->GetOuterRadiusUp();
292 sprintf(name,"US%2.2d",i);
294 sprintf(title,"TPC upper sector %d",i);
296 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
297 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
298 tubs->SetNumberOfDivisions(1);
300 nNode = new TNode(name,title,name,0,0,0,"");
301 nNode->SetLineColor(kColorTPC);
307 //_____________________________________________________________________________
308 void AliTPC::CreateMaterials()
310 //-----------------------------------------------
311 // Create Materials for for TPC simulations
312 //-----------------------------------------------
314 //-----------------------------------------------------------------
315 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
316 //-----------------------------------------------------------------
318 Int_t iSXFLD=gAlice->Field()->Integ();
319 Float_t sXMGMX=gAlice->Field()->Max();
321 Float_t amat[5]; // atomic numbers
322 Float_t zmat[5]; // z
323 Float_t wmat[5]; // proportions
329 //***************** Gases *************************
332 //--------------------------------------------------------------
333 // gases - air and CO2
334 //--------------------------------------------------------------
350 AliMixture(10,"CO2",amat,zmat,density,2,wmat);
365 AliMixture(11,"Air",amat,zmat,density,2,wmat);
367 //----------------------------------------------------------------
369 //----------------------------------------------------------------
372 // Drift gases 1 - nonsensitive, 2 - sensitive
373 // Ne-CO2-N (85-10-5)
392 AliMixture(12,"Ne-CO2-N-1",amat,zmat,density,4,wmat);
393 AliMixture(13,"Ne-CO2-N-2",amat,zmat,density,4,wmat);
394 AliMixture(30,"Ne-CO2-N-3",amat,zmat,density,4,wmat);
395 //----------------------------------------------------------------------
397 //----------------------------------------------------------------------
419 AliMixture(14,"Kevlar",amat,zmat,density,-4,wmat);
440 AliMixture(15,"NOMEX",amat,zmat,density,-4,wmat);
458 AliMixture(16,"Makrolon",amat,zmat,density,-3,wmat);
476 AliMixture(17, "Tedlar",amat,zmat,density,-3,wmat);
494 AliMixture(18, "Mylar",amat,zmat,density,-3,wmat);
495 // material for "prepregs"
496 // Epoxy - C14 H20 O3
499 // prepreg1 60% C-fiber, 40% epoxy (vol)
514 AliMixture(19, "Prepreg1",amat,zmat,density,3,wmat);
516 //prepreg2 60% glass-fiber, 40% epoxy
535 AliMixture(20, "Prepreg2",amat,zmat,density,4,wmat);
537 //prepreg3 50% glass-fiber, 50% epoxy
556 AliMixture(21, "Prepreg3",amat,zmat,density,4,wmat);
558 // G10 60% SiO2 40% epoxy
577 AliMixture(22, "G10",amat,zmat,density,4,wmat);
586 AliMaterial(23,"Al",amat[0],zmat[0],density,999.,999.);
588 // Si (for electronics
595 AliMaterial(24,"Si",amat[0],zmat[0],density,999.,999.);
604 AliMaterial(25,"Cu",amat[0],zmat[0],density,999.,999.);
606 // Epoxy - C14 H20 O3
622 AliMixture(26,"Epoxy",amat,zmat,density,-3,wmat);
640 AliMixture(27,"Plexiglas",amat,zmat,density,-3,wmat);
648 AliMaterial(28,"C",amat[0],zmat[0],density,999.,999.);
650 // Fe (steel for the inner heat screen)
658 AliMaterial(29,"Fe",amat[0],zmat[0],density,999.,999.);
660 // Peek - (C6H4-O-OC6H4-O-C6H4-CO)n
675 AliMixture(30,"Peek",amat,zmat,density,-3,wmat);
690 AliMixture(31,"Alumina",amat,zmat,density,-2,wmat);
709 AliMixture(32,"Water",amat,zmat,density,-2,wmat);
711 //----------------------------------------------------------
712 // tracking media for gases
713 //----------------------------------------------------------
715 AliMedium(0, "Air", 11, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
716 AliMedium(1, "Ne-CO2-N-1", 12, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
717 AliMedium(2, "Ne-CO2-N-2", 13, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
718 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
719 AliMedium(20, "Ne-CO2-N-3", 30, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
720 //-----------------------------------------------------------
721 // tracking media for solids
722 //-----------------------------------------------------------
724 AliMedium(4,"Al",23,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
725 AliMedium(5,"Kevlar",14,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
726 AliMedium(6,"Nomex",15,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
727 AliMedium(7,"Makrolon",16,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
728 AliMedium(8,"Mylar",18,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
729 AliMedium(9,"Tedlar",17,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
731 AliMedium(10,"Prepreg1",19,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
732 AliMedium(11,"Prepreg2",20,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
733 AliMedium(12,"Prepreg3",21,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
734 AliMedium(13,"Epoxy",26,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
736 AliMedium(14,"Cu",25,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
737 AliMedium(15,"Si",24,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
738 AliMedium(16,"G10",22,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
739 AliMedium(17,"Plexiglas",27,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
740 AliMedium(18,"Steel",29,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
741 AliMedium(19,"Peek",30,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
742 AliMedium(21,"Alumina",31,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
743 AliMedium(22,"Water",32,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
746 void AliTPC::GenerNoise(Int_t tablesize)
749 //Generate table with noise
756 if (fNoiseTable) delete[] fNoiseTable;
757 fNoiseTable = new Float_t[tablesize];
758 fNoiseDepth = tablesize;
759 fCurrentNoise =0; //!index of the noise in the noise table
761 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
762 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
765 Float_t AliTPC::GetNoise()
767 // get noise from table
768 // if ((fCurrentNoise%10)==0)
769 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
770 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
771 return fNoiseTable[fCurrentNoise++];
772 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
776 Bool_t AliTPC::IsSectorActive(Int_t sec) const
779 // check if the sector is active
780 if (!fActiveSectors) return kTRUE;
781 else return fActiveSectors[sec];
784 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
786 // activate interesting sectors
787 SetTreeAddress();//just for security
788 if (!fActiveSectors) fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
789 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
790 for (Int_t i=0;i<n;i++)
791 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
795 void AliTPC::SetActiveSectors(Int_t flag)
798 // activate sectors which were hitted by tracks
800 SetTreeAddress();//just for security
801 if (fHitType==0) return; // if Clones hit - not short volume ID information
802 if (!fActiveSectors) fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
804 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
807 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
811 AliFatal("Can not find TreeH in folder");
814 if (fHitType>1) branch = TreeH()->GetBranch("TPC2");
815 else branch = TreeH()->GetBranch("TPC");
816 Stat_t ntracks = TreeH()->GetEntries();
817 // loop over all hits
818 AliDebug(1,Form("Got %d tracks",ntracks));
820 for(Int_t track=0;track<ntracks;track++) {
823 if (fTrackHits && fHitType&4) {
824 TBranch * br1 = TreeH()->GetBranch("fVolumes");
825 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
826 br1->GetEvent(track);
827 br2->GetEvent(track);
828 Int_t *volumes = fTrackHits->GetVolumes();
829 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++) {
830 if (volumes[j]>-1 && volumes[j]<fTPCParam->GetNSector()) {
831 fActiveSectors[volumes[j]]=kTRUE;
834 AliError(Form("Volume %d -> sector number %d is outside (0..%d)",
837 fTPCParam->GetNSector()));
843 // if (fTrackHitsOld && fHitType&2) {
844 // TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
845 // br->GetEvent(track);
846 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
847 // for (UInt_t j=0;j<ar->GetSize();j++){
848 // fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
857 //_____________________________________________________________________________
858 void AliTPC::Digits2Raw()
860 // convert digits of the current event to raw data
862 static const Int_t kThreshold = 0;
864 fLoader->LoadDigits();
865 TTree* digits = fLoader->TreeD();
867 AliError("No digits tree");
872 AliSimDigits* digrow = &digarr;
873 digits->GetBranch("Segment")->SetAddress(&digrow);
875 const char* fileName = "AliTPCDDL.dat";
876 AliTPCBuffer* buffer = new AliTPCBuffer(fileName);
880 // 2: txt files with digits
881 //BE CAREFUL, verbose level 2 MUST be used only for debugging and
882 //it is highly suggested to use this mode only for debugging digits files
883 //reasonably small, because otherwise the size of the txt files can reach
884 //quickly several MB wasting time and disk space.
885 buffer->SetVerbose(0);
887 Int_t nEntries = Int_t(digits->GetEntries());
888 Int_t previousSector = -1;
890 for (Int_t i = 0; i < nEntries; i++) {
893 fTPCParam->AdjustSectorRow(digarr.GetID(), sector, row);
894 if(previousSector != sector) {
896 previousSector = sector;
899 if (sector < 36) { //inner sector [0;35]
901 //the whole row is written into the output file
902 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
903 sector, subSector, row);
905 //only the pads in the range [37;48] are written into the output file
906 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 1,
907 sector, subSector, row);
909 //only the pads outside the range [37;48] are written into the output file
910 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 2,
911 sector, subSector, row);
914 } else { //outer sector [36;71]
915 if (row == 54) subSector = 2;
916 if ((row != 27) && (row != 76)) {
917 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
918 sector, subSector, row);
919 } else if (row == 27) {
920 //only the pads outside the range [43;46] are written into the output file
921 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 2,
922 sector, subSector, row);
924 //only the pads in the range [43;46] are written into the output file
925 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 1,
926 sector, subSector, row);
927 } else if (row == 76) {
928 //only the pads outside the range [33;88] are written into the output file
929 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 2,
930 sector, subSector, row);
932 //only the pads in the range [33;88] are written into the output file
933 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 1,
934 sector, subSector, row);
940 fLoader->UnloadDigits();
942 AliTPCDDLRawData rawWriter;
943 rawWriter.SetVerbose(0);
945 rawWriter.RawData(fileName);
946 gSystem->Unlink(fileName);
951 //_____________________________________________________________________________
952 Bool_t AliTPC::Raw2SDigits(AliRawReader* rawReader){
953 // Converts the TPC raw data into summable digits
954 // The method is used for merging simulated and
956 if (fLoader->TreeS() == 0x0 ) {
957 fLoader->MakeTree("S");
960 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
962 //setup TPCDigitsArray
963 if(GetDigitsArray()) delete GetDigitsArray();
965 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
966 arr->SetClass("AliSimDigits");
967 arr->Setup(fTPCParam);
968 arr->MakeTree(fLoader->TreeS());
972 // set zero suppression to "0"
973 fTPCParam->SetZeroSup(0);
976 const Int_t kmaxTime = fTPCParam->GetMaxTBin();
977 const Int_t kNIS = fTPCParam->GetNInnerSector();
978 const Int_t kNOS = fTPCParam->GetNOuterSector();
979 const Int_t kNS = kNIS + kNOS;
981 Short_t** allBins = NULL; //array which contains the data for one sector
983 for(Int_t iSector = 0; iSector < kNS; iSector++) {
985 Int_t nRows = fTPCParam->GetNRow(iSector);
986 Int_t nDDLs = 0, indexDDL = 0;
987 if (iSector < kNIS) {
989 indexDDL = iSector * 2;
993 indexDDL = (iSector-kNIS) * 4 + kNIS * 2;
996 // Loas the raw data for corresponding DDLs
998 AliTPCRawStream input(rawReader);
999 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
1001 // Alocate and init the array with the sector data
1002 allBins = new Short_t*[nRows];
1003 for (Int_t iRow = 0; iRow < nRows; iRow++) {
1004 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
1005 Int_t maxBin = kmaxTime*maxPad;
1006 allBins[iRow] = new Short_t[maxBin];
1007 memset(allBins[iRow],0,sizeof(Short_t)*maxBin);
1010 // Begin loop over altro data
1011 while (input.Next()) {
1013 if (input.GetSector() != iSector)
1014 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",iSector,input.GetSector()));
1016 Int_t iRow = input.GetRow();
1017 if (iRow < 0 || iRow >= nRows)
1018 AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
1019 iRow, 0, nRows -1));
1020 Int_t iPad = input.GetPad();
1022 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
1024 if (iPad < 0 || iPad >= maxPad)
1025 AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
1026 iPad, 0, maxPad -1));
1028 Int_t iTimeBin = input.GetTime();
1029 if ( iTimeBin < 0 || iTimeBin >= kmaxTime)
1030 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
1031 iTimeBin, 0, kmaxTime -1));
1033 Int_t maxBin = kmaxTime*maxPad;
1035 if (((iPad*kmaxTime+iTimeBin) >= maxBin) ||
1036 ((iPad*kmaxTime+iTimeBin) < 0))
1037 AliFatal(Form("Index outside the allowed range"
1038 " Sector=%d Row=%d Pad=%d Timebin=%d"
1039 " (Max.index=%d)",iSector,iRow,iPad,iTimeBin,maxBin));
1041 allBins[iRow][iPad*kmaxTime+iTimeBin] = input.GetSignal();
1043 } // End loop over altro data
1045 // Now fill the digits array
1046 if (fDigitsArray->GetTree()==0) {
1047 AliFatal("Tree not set in fDigitsArray");
1050 for (Int_t iRow = 0; iRow < nRows; iRow++) {
1051 AliDigits * dig = fDigitsArray->CreateRow(iSector,iRow);
1053 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
1054 for(Int_t iPad = 0; iPad < maxPad; iPad++) {
1055 for(Int_t iTimeBin = 0; iTimeBin < kmaxTime; iTimeBin++) {
1056 Short_t q = allBins[iRow][iPad*kmaxTime + iTimeBin];
1057 if (q <= 0) continue;
1059 dig->SetDigitFast((Short_t)q,iTimeBin,iPad);
1062 fDigitsArray->StoreRow(iSector,iRow);
1063 Int_t ndig = dig->GetDigitSize();
1066 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1067 iSector,iRow,ndig));
1069 fDigitsArray->ClearRow(iSector,iRow);
1071 } // end of the sector digitization
1073 for (Int_t iRow = 0; iRow < nRows; iRow++)
1074 delete [] allBins[iRow];
1080 fLoader->WriteSDigits("OVERWRITE");
1082 if(GetDigitsArray()) delete GetDigitsArray();
1083 SetDigitsArray(0x0);
1088 //______________________________________________________________________
1089 AliDigitizer* AliTPC::CreateDigitizer(AliRunDigitizer* manager) const
1091 return new AliTPCDigitizer(manager);
1094 void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
1096 //create digits from summable digits
1097 GenerNoise(500000); //create teble with noise
1099 //conect tree with sSDigits
1100 TTree *t = fLoader->TreeS();
1103 fLoader->LoadSDigits("READ");
1104 t = fLoader->TreeS();
1106 AliError("Can not get input TreeS");
1111 if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
1113 AliSimDigits digarr, *dummy=&digarr;
1114 TBranch* sdb = t->GetBranch("Segment");
1116 AliError("Can not find branch with segments in TreeS.");
1120 sdb->SetAddress(&dummy);
1122 Stat_t nentries = t->GetEntries();
1124 // set zero suppression
1126 fTPCParam->SetZeroSup(2);
1128 // get zero suppression
1130 Int_t zerosup = fTPCParam->GetZeroSup();
1132 //make tree with digits
1134 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1135 arr->SetClass("AliSimDigits");
1136 arr->Setup(fTPCParam);
1137 arr->MakeTree(fLoader->TreeD());
1139 AliTPCParam * par = fTPCParam;
1141 //Loop over segments of the TPC
1143 for (Int_t n=0; n<nentries; n++) {
1146 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1147 AliWarning(Form("Invalid segment ID ! %d",digarr.GetID()));
1150 if (!IsSectorActive(sec)) continue;
1152 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1153 Int_t nrows = digrow->GetNRows();
1154 Int_t ncols = digrow->GetNCols();
1156 digrow->ExpandBuffer();
1157 digarr.ExpandBuffer();
1158 digrow->ExpandTrackBuffer();
1159 digarr.ExpandTrackBuffer();
1162 Short_t * pamp0 = digarr.GetDigits();
1163 Int_t * ptracks0 = digarr.GetTracks();
1164 Short_t * pamp1 = digrow->GetDigits();
1165 Int_t * ptracks1 = digrow->GetTracks();
1166 Int_t nelems =nrows*ncols;
1167 Int_t saturation = fTPCParam->GetADCSat() - 1;
1168 //use internal structure of the AliDigits - for speed reason
1169 //if you cahnge implementation
1170 //of the Alidigits - it must be rewriten -
1171 for (Int_t i= 0; i<nelems; i++){
1172 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1174 if (q>saturation) q=saturation;
1177 ptracks1[0]=ptracks0[0];
1178 ptracks1[nelems]=ptracks0[nelems];
1179 ptracks1[2*nelems]=ptracks0[2*nelems];
1187 arr->StoreRow(sec,row);
1188 arr->ClearRow(sec,row);
1193 fLoader->WriteDigits("OVERWRITE");
1197 //__________________________________________________________________
1198 void AliTPC::SetDefaults(){
1200 // setting the defaults
1203 // Set response functions
1206 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1208 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1210 AliInfo("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...");
1212 param = new AliTPCParamSR();
1215 param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1218 AliFatal("No TPC parameters found");
1222 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1223 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1224 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1225 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1226 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1227 rf->SetOffset(3*param->GetZSigma());
1230 TDirectory *savedir=gDirectory;
1231 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1233 AliFatal("Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !");
1236 prfinner->Read("prf_07504_Gati_056068_d02");
1237 //PH Set different names
1238 s=prfinner->GetGRF()->GetName();
1240 prfinner->GetGRF()->SetName(s.Data());
1242 prfouter1->Read("prf_10006_Gati_047051_d03");
1243 s=prfouter1->GetGRF()->GetName();
1245 prfouter1->GetGRF()->SetName(s.Data());
1247 prfouter2->Read("prf_15006_Gati_047051_d03");
1248 s=prfouter2->GetGRF()->GetName();
1250 prfouter2->GetGRF()->SetName(s.Data());
1255 param->SetInnerPRF(prfinner);
1256 param->SetOuter1PRF(prfouter1);
1257 param->SetOuter2PRF(prfouter2);
1258 param->SetTimeRF(rf);
1268 //__________________________________________________________________
1269 void AliTPC::Hits2Digits()
1272 // creates digits from hits
1274 if (!fTPCParam->IsGeoRead()){
1276 // read transformation matrices for gGeoManager
1278 fTPCParam->ReadGeoMatrices();
1281 fLoader->LoadHits("read");
1282 fLoader->LoadDigits("recreate");
1283 AliRunLoader* runLoader = fLoader->GetRunLoader();
1285 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1286 //PH runLoader->GetEvent(iEvent);
1287 Hits2Digits(iEvent);
1290 fLoader->UnloadHits();
1291 fLoader->UnloadDigits();
1293 //__________________________________________________________________
1294 void AliTPC::Hits2Digits(Int_t eventnumber)
1296 //----------------------------------------------------
1297 // Loop over all sectors for a single event
1298 //----------------------------------------------------
1299 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1300 rl->GetEvent(eventnumber);
1302 if (fLoader->TreeH() == 0x0) {
1303 if(fLoader->LoadHits()) {
1304 AliError("Can not load hits.");
1309 if (fLoader->TreeD() == 0x0 ) {
1310 fLoader->MakeTree("D");
1311 if (fLoader->TreeD() == 0x0 ) {
1312 AliError("Can not get TreeD");
1317 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1318 GenerNoise(500000); //create teble with noise
1320 //setup TPCDigitsArray
1322 if(GetDigitsArray()) delete GetDigitsArray();
1324 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1325 arr->SetClass("AliSimDigits");
1326 arr->Setup(fTPCParam);
1328 arr->MakeTree(fLoader->TreeD());
1329 SetDigitsArray(arr);
1331 fDigitsSwitch=0; // standard digits
1333 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1334 if (IsSectorActive(isec)) {
1335 AliDebug(1,Form("Hits2Digits","Sector %d is active.",isec));
1336 Hits2DigitsSector(isec);
1339 AliDebug(1,Form("Hits2Digits","Sector %d is NOT active.",isec));
1342 fLoader->WriteDigits("OVERWRITE");
1344 //this line prevents the crash in the similar one
1345 //on the beginning of this method
1346 //destructor attempts to reset the tree, which is deleted by the loader
1347 //need to be redesign
1348 if(GetDigitsArray()) delete GetDigitsArray();
1349 SetDigitsArray(0x0);
1353 //__________________________________________________________________
1354 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1357 //-----------------------------------------------------------
1358 // summable digits - 16 bit "ADC", no noise, no saturation
1359 //-----------------------------------------------------------
1361 //----------------------------------------------------
1362 // Loop over all sectors for a single event
1363 //----------------------------------------------------
1365 AliRunLoader* rl = fLoader->GetRunLoader();
1367 rl->GetEvent(eventnumber);
1368 if (fLoader->TreeH() == 0x0) {
1369 if(fLoader->LoadHits()) {
1370 AliError("Can not load hits.");
1377 if (fLoader->TreeS() == 0x0 ) {
1378 fLoader->MakeTree("S");
1381 if(fDefaults == 0) SetDefaults();
1383 GenerNoise(500000); //create table with noise
1384 //setup TPCDigitsArray
1386 if(GetDigitsArray()) delete GetDigitsArray();
1389 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1390 arr->SetClass("AliSimDigits");
1391 arr->Setup(fTPCParam);
1392 arr->MakeTree(fLoader->TreeS());
1394 SetDigitsArray(arr);
1396 fDigitsSwitch=1; // summable digits
1398 // set zero suppression to "0"
1400 fTPCParam->SetZeroSup(0);
1402 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1403 if (IsSectorActive(isec)) {
1404 Hits2DigitsSector(isec);
1407 fLoader->WriteSDigits("OVERWRITE");
1409 //this line prevents the crash in the similar one
1410 //on the beginning of this method
1411 //destructor attempts to reset the tree, which is deleted by the loader
1412 //need to be redesign
1413 if(GetDigitsArray()) delete GetDigitsArray();
1414 SetDigitsArray(0x0);
1416 //__________________________________________________________________
1418 void AliTPC::Hits2SDigits()
1421 //-----------------------------------------------------------
1422 // summable digits - 16 bit "ADC", no noise, no saturation
1423 //-----------------------------------------------------------
1425 if (!fTPCParam->IsGeoRead()){
1427 // read transformation matrices for gGeoManager
1429 fTPCParam->ReadGeoMatrices();
1432 fLoader->LoadHits("read");
1433 fLoader->LoadSDigits("recreate");
1434 AliRunLoader* runLoader = fLoader->GetRunLoader();
1436 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1437 runLoader->GetEvent(iEvent);
1440 Hits2SDigits2(iEvent);
1443 fLoader->UnloadHits();
1444 fLoader->UnloadSDigits();
1446 //_____________________________________________________________________________
1448 void AliTPC::Hits2DigitsSector(Int_t isec)
1450 //-------------------------------------------------------------------
1451 // TPC conversion from hits to digits.
1452 //-------------------------------------------------------------------
1454 //-----------------------------------------------------------------
1455 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1456 //-----------------------------------------------------------------
1458 //-------------------------------------------------------
1459 // Get the access to the track hits
1460 //-------------------------------------------------------
1462 // check if the parameters are set - important if one calls this method
1463 // directly, not from the Hits2Digits
1465 if(fDefaults == 0) SetDefaults();
1467 TTree *tH = TreeH(); // pointer to the hits tree
1469 AliFatal("Can not find TreeH in folder");
1473 Stat_t ntracks = tH->GetEntries();
1477 //-------------------------------------------
1478 // Only if there are any tracks...
1479 //-------------------------------------------
1483 Int_t nrows =fTPCParam->GetNRow(isec);
1485 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1487 MakeSector(isec,nrows,tH,ntracks,row);
1489 //--------------------------------------------------------
1490 // Digitize this sector, row by row
1491 // row[i] is the pointer to the TObjArray of TVectors,
1492 // each one containing electrons accepted on this
1493 // row, assigned into tracks
1494 //--------------------------------------------------------
1498 if (fDigitsArray->GetTree()==0) {
1499 AliFatal("Tree not set in fDigitsArray");
1502 for (i=0;i<nrows;i++){
1504 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1506 DigitizeRow(i,isec,row);
1508 fDigitsArray->StoreRow(isec,i);
1510 Int_t ndig = dig->GetDigitSize();
1513 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1516 fDigitsArray->ClearRow(isec,i);
1519 } // end of the sector digitization
1521 for(i=0;i<nrows+2;i++){
1526 delete [] row; // delete the array of pointers to TObjArray-s
1530 } // end of Hits2DigitsSector
1533 //_____________________________________________________________________________
1534 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1536 //-----------------------------------------------------------
1537 // Single row digitization, coupling from the neighbouring
1538 // rows taken into account
1539 //-----------------------------------------------------------
1541 //-----------------------------------------------------------------
1542 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1543 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1544 //-----------------------------------------------------------------
1546 Float_t zerosup = fTPCParam->GetZeroSup();
1548 fCurrentIndex[1]= isec;
1551 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1552 Int_t nofTbins = fTPCParam->GetMaxTBin();
1553 Int_t indexRange[4];
1555 // Integrated signal for this row
1556 // and a single track signal
1559 TMatrixF *m1 = new TMatrixF(0,nofPads,0,nofTbins); // integrated
1560 TMatrixF *m2 = new TMatrixF(0,nofPads,0,nofTbins); // single
1562 TMatrixF &total = *m1;
1564 // Array of pointers to the label-signal list
1566 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1567 Float_t **pList = new Float_t* [nofDigits];
1571 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1577 for (Int_t row= row1;row<=row2;row++){
1578 Int_t nTracks= rows[row]->GetEntries();
1579 for (i1=0;i1<nTracks;i1++){
1580 fCurrentIndex[2]= row;
1581 fCurrentIndex[3]=irow+1;
1583 m2->Zero(); // clear single track signal matrix
1584 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1585 GetList(trackLabel,nofPads,m2,indexRange,pList);
1587 else GetSignal(rows[row],i1,0,m1,indexRange);
1593 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1595 Float_t fzerosup = zerosup+0.5;
1596 for(Int_t it=0;it<nofTbins;it++){
1597 for(Int_t ip=0;ip<nofPads;ip++){
1599 Float_t q=total(ip,it);
1600 if(fDigitsSwitch == 0){
1602 if(q <=fzerosup) continue; // do not fill zeros
1604 if(q >= fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat() - 1; // saturation
1609 if(q <= 0.) continue; // do not fill zeros
1610 if(q>2000.) q=2000.;
1616 // "real" signal or electronic noise (list = -1)?
1619 for(Int_t j1=0;j1<3;j1++){
1620 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1625 <A NAME="AliDigits"></A>
1626 using of AliDigits object
1629 dig->SetDigitFast((Short_t)q,it,ip);
1630 if (fDigitsArray->IsSimulated()) {
1631 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1632 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1633 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1636 } // end of loop over time buckets
1637 } // end of lop over pads
1640 // This row has been digitized, delete nonused stuff
1643 for(lp=0;lp<nofDigits;lp++){
1644 if(pList[lp]) delete [] pList[lp];
1652 } // end of DigitizeRow
1654 //_____________________________________________________________________________
1656 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1657 TMatrixF *m1, TMatrixF *m2,Int_t *indexRange)
1660 //---------------------------------------------------------------
1661 // Calculates 2-D signal (pad,time) for a single track,
1662 // returns a pointer to the signal matrix and the track label
1663 // No digitization is performed at this level!!!
1664 //---------------------------------------------------------------
1666 //-----------------------------------------------------------------
1667 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1668 // Modified: Marian Ivanov
1669 //-----------------------------------------------------------------
1673 tv = (TVector*)p1->At(ntr); // pointer to a track
1676 Float_t label = v(0);
1677 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1))/2;
1679 Int_t nElectrons = (tv->GetNrows()-1)/5;
1680 indexRange[0]=9999; // min pad
1681 indexRange[1]=-1; // max pad
1682 indexRange[2]=9999; //min time
1683 indexRange[3]=-1; // max time
1685 TMatrixF &signal = *m1;
1686 TMatrixF &total = *m2;
1688 // Loop over all electrons
1690 for(Int_t nel=0; nel<nElectrons; nel++){
1692 Float_t aval = v(idx+4);
1693 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1694 Float_t xyz[4]={v(idx+1),v(idx+2),v(idx+3),v(idx+5)};
1695 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
1697 Int_t *index = fTPCParam->GetResBin(0);
1698 Float_t *weight = & (fTPCParam->GetResWeight(0));
1700 if (n>0) for (Int_t i =0; i<n; i++){
1701 //Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1702 Int_t pad = (xyz[1]>0.) ? index[1]+centralPad : index[1]+centralPad -1; //corr M.K.
1703 pad = (xyz[2]>0.) ? pad : 2*centralPad-pad-1; // corr M.K.
1707 Int_t time=index[2];
1708 Float_t qweight = *(weight)*eltoadcfac;
1710 if (m1!=0) signal(pad,time)+=qweight;
1711 total(pad,time)+=qweight;
1712 if (indexRange[0]>pad) indexRange[0]=pad;
1713 if (indexRange[1]<pad) indexRange[1]=pad;
1714 if (indexRange[2]>time) indexRange[2]=time;
1715 if (indexRange[3]<time) indexRange[3]=time;
1722 } // end of loop over electrons
1724 return label; // returns track label when finished
1727 //_____________________________________________________________________________
1728 void AliTPC::GetList(Float_t label,Int_t np,TMatrixF *m,
1729 Int_t *indexRange, Float_t **pList)
1731 //----------------------------------------------------------------------
1732 // Updates the list of tracks contributing to digits for a given row
1733 //----------------------------------------------------------------------
1735 //-----------------------------------------------------------------
1736 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1737 //-----------------------------------------------------------------
1739 TMatrixF &signal = *m;
1741 // lop over nonzero digits
1743 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1744 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1747 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1749 if(signal(ip,it)<0.5) continue;
1751 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1753 if(!pList[globalIndex]){
1756 // Create new list (6 elements - 3 signals and 3 labels),
1759 pList[globalIndex] = new Float_t [6];
1763 *pList[globalIndex] = -1.;
1764 *(pList[globalIndex]+1) = -1.;
1765 *(pList[globalIndex]+2) = -1.;
1766 *(pList[globalIndex]+3) = -1.;
1767 *(pList[globalIndex]+4) = -1.;
1768 *(pList[globalIndex]+5) = -1.;
1770 *pList[globalIndex] = label;
1771 *(pList[globalIndex]+3) = signal(ip,it);
1775 // check the signal magnitude
1777 Float_t highest = *(pList[globalIndex]+3);
1778 Float_t middle = *(pList[globalIndex]+4);
1779 Float_t lowest = *(pList[globalIndex]+5);
1782 // compare the new signal with already existing list
1785 if(signal(ip,it)<lowest) continue; // neglect this track
1789 if (signal(ip,it)>highest){
1790 *(pList[globalIndex]+5) = middle;
1791 *(pList[globalIndex]+4) = highest;
1792 *(pList[globalIndex]+3) = signal(ip,it);
1794 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1795 *(pList[globalIndex]+1) = *pList[globalIndex];
1796 *pList[globalIndex] = label;
1798 else if (signal(ip,it)>middle){
1799 *(pList[globalIndex]+5) = middle;
1800 *(pList[globalIndex]+4) = signal(ip,it);
1802 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1803 *(pList[globalIndex]+1) = label;
1806 *(pList[globalIndex]+5) = signal(ip,it);
1807 *(pList[globalIndex]+2) = label;
1811 } // end of loop over pads
1812 } // end of loop over time bins
1815 //___________________________________________________________________
1816 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1817 Stat_t ntracks,TObjArray **row)
1820 //-----------------------------------------------------------------
1821 // Prepares the sector digitization, creates the vectors of
1822 // tracks for each row of this sector. The track vector
1823 // contains the track label and the position of electrons.
1824 //-----------------------------------------------------------------
1827 // The trasport of the electrons through TPC drift volume
1828 // Drift (drift velocity + velocity map(not yet implemented)))
1829 // Application of the random processes (diffusion, gas gain)
1830 // Systematic effects (ExB effect in drift volume + ROCs)
1833 // Loop over primary electrons:
1834 // Creation of the secondary electrons
1835 // Loop over electrons (primary+ secondaries)
1836 // Global coordinate frame:
1837 // 1. Skip electrons if attached
1838 // 2. ExB effect in drift volume
1839 // a.) Simulation calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1840 // b.) Reconstruction - calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1841 // 3. Generation of gas gain (Random - Exponential distribution)
1842 // 4. TransportElectron function (diffusion)
1844 // 5. Conversion to the local coordinate frame pad-row, pad, timebin
1845 // 6. Apply Time0 shift - AliTPCCalPad class
1846 // a.) Plus sign in simulation
1847 // b.) Minus sign in reconstruction
1850 //-----------------------------------------------------------------
1851 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1852 // Origin: Marian Ivanov, marian.ivanov@cern.ch
1853 //-----------------------------------------------------------------
1854 AliTPCcalibDB* const calib=AliTPCcalibDB::Instance();
1856 Float_t gasgain = fTPCParam->GetGasGain();
1857 gasgain = gasgain/fGainFactor;
1861 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1864 if (fHitType>1) branch = TH->GetBranch("TPC2");
1865 else branch = TH->GetBranch("TPC");
1868 //----------------------------------------------
1869 // Create TObjArray-s, one for each row,
1870 // each TObjArray will store the TVectors
1871 // of electrons, one TVectors per each track.
1872 //----------------------------------------------
1874 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
1875 TVector **tracks = new TVector* [nrows+2]; //pointers to the track vectors
1877 for(i=0; i<nrows+2; i++){
1878 row[i] = new TObjArray;
1885 //--------------------------------------------------------------------
1886 // Loop over tracks, the "track" contains the full history
1887 //--------------------------------------------------------------------
1889 Int_t previousTrack,currentTrack;
1890 previousTrack = -1; // nothing to store so far!
1892 for(Int_t track=0;track<ntracks;track++){
1893 Bool_t isInSector=kTRUE;
1895 isInSector = TrackInVolume(isec,track);
1896 if (!isInSector) continue;
1898 branch->GetEntry(track); // get next track
1902 tpcHit = (AliTPChit*)FirstHit(-1);
1904 //--------------------------------------------------------------
1906 //--------------------------------------------------------------
1911 Int_t sector=tpcHit->fSector; // sector number
1913 tpcHit = (AliTPChit*) NextHit();
1917 // Remove hits which arrive before the TPC opening gate signal
1918 if(((fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
1919 /fTPCParam->GetDriftV()+tpcHit->Time())<fTPCParam->GetGateDelay()) {
1920 tpcHit = (AliTPChit*) NextHit();
1924 currentTrack = tpcHit->Track(); // track number
1926 if(currentTrack != previousTrack){
1928 // store already filled fTrack
1930 for(i=0;i<nrows+2;i++){
1931 if(previousTrack != -1){
1932 if(nofElectrons[i]>0){
1933 TVector &v = *tracks[i];
1934 v(0) = previousTrack;
1935 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
1936 row[i]->Add(tracks[i]);
1939 delete tracks[i]; // delete empty TVector
1945 tracks[i] = new TVector(601); // TVectors for the next fTrack
1947 } // end of loop over rows
1949 previousTrack=currentTrack; // update track label
1952 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1954 //---------------------------------------------------
1955 // Calculate the electron attachment probability
1956 //---------------------------------------------------
1959 Float_t time = 1.e6*(fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
1960 /fTPCParam->GetDriftV();
1962 Float_t attProb = fTPCParam->GetAttCoef()*
1963 fTPCParam->GetOxyCont()*time; // fraction!
1965 //-----------------------------------------------
1966 // Loop over electrons
1967 //-----------------------------------------------
1970 for(Int_t nel=0;nel<qI;nel++){
1971 // skip if electron lost due to the attachment
1972 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1977 Double_t dxyz0[3],dxyz1[3];
1978 dxyz0[0]=tpcHit->X();
1979 dxyz0[1]=tpcHit->Y();
1980 dxyz0[2]=tpcHit->Z();
1981 if (calib->GetExB()){
1982 calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1984 AliError("Not valid ExB calibration");
1985 dxyz1[0]=tpcHit->X();
1986 dxyz1[1]=tpcHit->Y();
1987 dxyz1[2]=tpcHit->Z();
1995 // protection for the nonphysical avalanche size (10**6 maximum)
1997 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
1998 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2001 TransportElectron(xyz,index);
2003 fTPCParam->GetPadRow(xyz,index);
2005 // Add Time0 correction due unisochronity
2006 // xyz[0] - pad row coordinate
2007 // xyz[1] - pad coordinate
2008 // xyz[2] - is in now time bin coordinate system
2009 Float_t correction =0;
2010 if (calib->GetPadTime0()){
2011 if (!calib->GetPadTime0()->GetCalROC(isec)) continue;
2012 correction = calib->GetPadTime0()->GetCalROC(isec)->GetValue(TMath::Nint(xyz[0]),TMath::Nint(xyz[1]));
2016 // Electron track time (for pileup simulation)
2017 xyz[4] = tpcHit->Time()/fTPCParam->GetTSample();
2018 // row 0 - cross talk from the innermost row
2019 // row fNRow+1 cross talk from the outermost row
2020 rowNumber = index[2]+1;
2021 //transform position to local digit coordinates
2022 //relative to nearest pad row
2023 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2025 if (isec <fTPCParam->GetNInnerSector()) {
2026 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2027 y1 = fTPCParam->GetYInner(rowNumber);
2030 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2031 y1 = fTPCParam->GetYOuter(rowNumber);
2033 // gain inefficiency at the wires edges - linear
2036 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));
2038 nofElectrons[rowNumber]++;
2039 //----------------------------------
2040 // Expand vector if necessary
2041 //----------------------------------
2042 if(nofElectrons[rowNumber]>120){
2043 Int_t range = tracks[rowNumber]->GetNrows();
2044 if((nofElectrons[rowNumber])>(range-1)/5){
2046 tracks[rowNumber]->ResizeTo(range+500); // Add 100 electrons
2050 TVector &v = *tracks[rowNumber];
2051 Int_t idx = 5*nofElectrons[rowNumber]-4;
2052 Real_t * position = &(((TVector&)v)(idx)); //make code faster
2053 memcpy(position,xyz,5*sizeof(Float_t));
2055 } // end of loop over electrons
2057 tpcHit = (AliTPChit*)NextHit();
2059 } // end of loop over hits
2060 } // end of loop over tracks
2063 // store remaining track (the last one) if not empty
2066 for(i=0;i<nrows+2;i++){
2067 if(nofElectrons[i]>0){
2068 TVector &v = *tracks[i];
2069 v(0) = previousTrack;
2070 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
2071 row[i]->Add(tracks[i]);
2080 delete [] nofElectrons;
2082 } // end of MakeSector
2085 //_____________________________________________________________________________
2089 // Initialise TPC detector after definition of geometry
2091 AliDebug(1,"*********************************************");
2094 //_____________________________________________________________________________
2095 void AliTPC::MakeBranch(Option_t* option)
2098 // Create Tree branches for the TPC.
2101 Int_t buffersize = 4000;
2102 char branchname[10];
2103 sprintf(branchname,"%s",GetName());
2105 const char *h = strstr(option,"H");
2107 if ( h && (fHitType<=1) && (fHits == 0x0)) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2109 AliDetector::MakeBranch(option);
2111 const char *d = strstr(option,"D");
2113 if (fDigits && fLoader->TreeD() && d) {
2114 MakeBranchInTree(gAlice->TreeD(), branchname, &fDigits, buffersize, 0);
2117 if (fHitType>1) MakeBranch2(option,0); // MI change 14.09.2000
2120 //_____________________________________________________________________________
2121 void AliTPC::ResetDigits()
2124 // Reset number of digits and the digits array for this detector
2127 if (fDigits) fDigits->Clear();
2132 //_____________________________________________________________________________
2133 void AliTPC::SetSens(Int_t sens)
2136 //-------------------------------------------------------------
2137 // Activates/deactivates the sensitive strips at the center of
2138 // the pad row -- this is for the space-point resolution calculations
2139 //-------------------------------------------------------------
2141 //-----------------------------------------------------------------
2142 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2143 //-----------------------------------------------------------------
2149 void AliTPC::SetSide(Float_t side=0.)
2151 // choice of the TPC side
2156 //_____________________________________________________________________________
2158 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2161 // electron transport taking into account:
2163 // 2.ExB at the wires
2164 // 3. nonisochronity
2166 // xyz and index must be already transformed to system 1
2169 fTPCParam->Transform1to2(xyz,index);
2172 Float_t driftl=xyz[2];
2173 if(driftl<0.01) driftl=0.01;
2174 driftl=TMath::Sqrt(driftl);
2175 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2176 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2177 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2178 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2179 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2183 if (fTPCParam->GetMWPCReadout()==kTRUE){
2184 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2185 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2187 //add nonisochronity (not implemented yet)
2193 //______________________________________________________________________
2194 AliTPChit::AliTPChit()
2206 //_____________________________________________________________________________
2207 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits)
2208 :AliHit(shunt,track),
2215 // Creates a TPC hit object
2226 //________________________________________________________________________
2227 // Additional code because of the AliTPCTrackHitsV2
2229 void AliTPC::MakeBranch2(Option_t *option,const char */*file*/)
2232 // Create a new branch in the current Root Tree
2233 // The branch of fHits is automatically split
2234 // MI change 14.09.2000
2236 if (fHitType<2) return;
2237 char branchname[10];
2238 sprintf(branchname,"%s2",GetName());
2240 // Get the pointer to the header
2241 const char *cH = strstr(option,"H");
2243 if (fTrackHits && TreeH() && cH && fHitType&4) {
2244 AliDebug(1,"Making branch for Type 4 Hits");
2245 TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2248 // if (fTrackHitsOld && TreeH() && cH && fHitType&2) {
2249 // AliDebug(1,"Making branch for Type 2 Hits");
2250 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2251 // TreeH(),fBufferSize,99);
2252 // TreeH()->GetListOfBranches()->Add(branch);
2256 void AliTPC::SetTreeAddress()
2258 //Sets tree address for hits
2260 if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2261 AliDetector::SetTreeAddress();
2263 if (fHitType>1) SetTreeAddress2();
2266 void AliTPC::SetTreeAddress2()
2269 // Set branch address for the TrackHits Tree
2274 char branchname[20];
2275 sprintf(branchname,"%s2",GetName());
2277 // Branch address for hit tree
2278 TTree *treeH = TreeH();
2279 if ((treeH)&&(fHitType&4)) {
2280 branch = treeH->GetBranch(branchname);
2282 branch->SetAddress(&fTrackHits);
2283 AliDebug(1,"fHitType&4 Setting");
2286 AliDebug(1,"fHitType&4 Failed (can not find branch)");
2289 // if ((treeH)&&(fHitType&2)) {
2290 // branch = treeH->GetBranch(branchname);
2292 // branch->SetAddress(&fTrackHitsOld);
2293 // AliDebug(1,"fHitType&2 Setting");
2296 // AliDebug(1,"fHitType&2 Failed (can not find branch)");
2300 void AliTPC::FinishPrimary()
2302 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2303 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2307 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2310 // add hit to the list
2313 int primary = gAlice->GetMCApp()->GetPrimary(track);
2314 gAlice->GetMCApp()->Particle(primary)->SetBit(kKeepBit);
2318 gAlice->GetMCApp()->FlagTrack(track);
2320 if (fTrackHits && fHitType&4)
2321 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2322 hits[1],hits[2],(Int_t)hits[3],hits[4]);
2323 // if (fTrackHitsOld &&fHitType&2 )
2324 // fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2325 // hits[1],hits[2],(Int_t)hits[3]);
2329 void AliTPC::ResetHits()
2331 if (fHitType&1) AliDetector::ResetHits();
2332 if (fHitType>1) ResetHits2();
2335 void AliTPC::ResetHits2()
2339 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2340 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2344 AliHit* AliTPC::FirstHit(Int_t track)
2346 if (fHitType>1) return FirstHit2(track);
2347 return AliDetector::FirstHit(track);
2349 AliHit* AliTPC::NextHit()
2354 if (fHitType>1) return NextHit2();
2356 return AliDetector::NextHit();
2359 AliHit* AliTPC::FirstHit2(Int_t track)
2362 // Initialise the hit iterator
2363 // Return the address of the first hit for track
2364 // If track>=0 the track is read from disk
2365 // while if track<0 the first hit of the current
2366 // track is returned
2369 gAlice->ResetHits();
2370 TreeH()->GetEvent(track);
2373 if (fTrackHits && fHitType&4) {
2374 fTrackHits->First();
2375 return fTrackHits->GetHit();
2377 // if (fTrackHitsOld && fHitType&2) {
2378 // fTrackHitsOld->First();
2379 // return fTrackHitsOld->GetHit();
2385 AliHit* AliTPC::NextHit2()
2388 //Return the next hit for the current track
2391 // if (fTrackHitsOld && fHitType&2) {
2392 // fTrackHitsOld->Next();
2393 // return fTrackHitsOld->GetHit();
2397 return fTrackHits->GetHit();
2403 void AliTPC::LoadPoints(Int_t)
2408 if(fHitType==1) AliDetector::LoadPoints(a);
2409 else LoadPoints2(a);
2413 void AliTPC::RemapTrackHitIDs(Int_t *map)
2418 if (!fTrackHits) return;
2420 // if (fTrackHitsOld && fHitType&2){
2421 // AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2422 // for (UInt_t i=0;i<arr->GetSize();i++){
2423 // AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2424 // info->fTrackID = map[info->fTrackID];
2427 // if (fTrackHitsOld && fHitType&4){
2428 if (fTrackHits && fHitType&4){
2429 TClonesArray * arr = fTrackHits->GetArray();;
2430 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2431 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2432 info->SetTrackID(map[info->GetTrackID()]);
2437 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2439 //return bool information - is track in given volume
2440 //load only part of the track information
2441 //return true if current track is in volume
2444 // if (fTrackHitsOld && fHitType&2) {
2445 // TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
2446 // br->GetEvent(track);
2447 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2448 // for (UInt_t j=0;j<ar->GetSize();j++){
2449 // if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2453 if (fTrackHits && fHitType&4) {
2454 TBranch * br1 = TreeH()->GetBranch("fVolumes");
2455 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
2456 br2->GetEvent(track);
2457 br1->GetEvent(track);
2458 Int_t *volumes = fTrackHits->GetVolumes();
2459 Int_t nvolumes = fTrackHits->GetNVolumes();
2460 if (!volumes && nvolumes>0) {
2461 AliWarning(Form("Problematic track\t%d\t%d",track,nvolumes));
2464 for (Int_t j=0;j<nvolumes; j++)
2465 if (volumes[j]==id) return kTRUE;
2469 TBranch * br = TreeH()->GetBranch("fSector");
2470 br->GetEvent(track);
2471 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2472 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2479 //_____________________________________________________________________________
2480 void AliTPC::LoadPoints2(Int_t)
2483 // Store x, y, z of all hits in memory
2485 // if (fTrackHits == 0 && fTrackHitsOld==0) return;
2486 if (fTrackHits == 0 ) return;
2489 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2490 // if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2492 if (nhits == 0) return;
2493 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2494 if (fPoints == 0) fPoints = new TObjArray(tracks);
2497 Int_t *ntrk=new Int_t[tracks];
2498 Int_t *limi=new Int_t[tracks];
2499 Float_t **coor=new Float_t*[tracks];
2500 for(Int_t i=0;i<tracks;i++) {
2506 AliPoints *points = 0;
2509 Int_t chunk=nhits/4+1;
2511 // Loop over all the hits and store their position
2513 ahit = FirstHit2(-1);
2515 trk=ahit->GetTrack();
2516 if(ntrk[trk]==limi[trk]) {
2518 // Initialise a new track
2519 fp=new Float_t[3*(limi[trk]+chunk)];
2521 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2522 delete [] coor[trk];
2529 fp[3*ntrk[trk] ] = ahit->X();
2530 fp[3*ntrk[trk]+1] = ahit->Y();
2531 fp[3*ntrk[trk]+2] = ahit->Z();
2539 for(trk=0; trk<tracks; ++trk) {
2541 points = new AliPoints();
2542 points->SetMarkerColor(kYellow); //PH kYellow it the default in TPC
2543 points->SetMarkerSize(1);//PH Default size=1
2544 points->SetDetector(this);
2545 points->SetParticle(trk);
2546 points->SetPolyMarker(ntrk[trk],coor[trk],1); // Default style=1
2547 fPoints->AddAt(points,trk);
2548 delete [] coor[trk];
2558 //_____________________________________________________________________________
2559 void AliTPC::LoadPoints3(Int_t)
2562 // Store x, y, z of all hits in memory
2563 // - only intersection point with pad row
2564 if (fTrackHits == 0) return;
2566 Int_t nhits = fTrackHits->GetEntriesFast();
2567 if (nhits == 0) return;
2568 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2569 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2570 fPoints->Expand(2*tracks);
2573 Int_t *ntrk=new Int_t[tracks];
2574 Int_t *limi=new Int_t[tracks];
2575 Float_t **coor=new Float_t*[tracks];
2576 for(Int_t i=0;i<tracks;i++) {
2582 AliPoints *points = 0;
2585 Int_t chunk=nhits/4+1;
2587 // Loop over all the hits and store their position
2589 ahit = FirstHit2(-1);
2593 trk=ahit->GetTrack();
2594 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2595 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2596 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2597 if (currentrow!=lastrow){
2598 lastrow = currentrow;
2599 //later calculate intersection point
2600 if(ntrk[trk]==limi[trk]) {
2602 // Initialise a new track
2603 fp=new Float_t[3*(limi[trk]+chunk)];
2605 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2606 delete [] coor[trk];
2613 fp[3*ntrk[trk] ] = ahit->X();
2614 fp[3*ntrk[trk]+1] = ahit->Y();
2615 fp[3*ntrk[trk]+2] = ahit->Z();
2622 for(trk=0; trk<tracks; ++trk) {
2624 points = new AliPoints();
2625 points->SetMarkerColor(kMagenta); //PH kYellow + 1 is kMagenta
2626 points->SetMarkerStyle(5);
2627 points->SetMarkerSize(0.2);
2628 points->SetDetector(this);
2629 points->SetParticle(trk);
2630 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2631 fPoints->AddAt(points,tracks+trk);
2632 delete [] coor[trk];
2643 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
2646 fLoader = new AliTPCLoader(GetName(),topfoldername);
2650 ////////////////////////////////////////////////////////////////////////
2651 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
2653 // load TPC paarmeters from a given file or create new if the object
2654 // is not found there
2655 // 12/05/2003 This method should be moved to the AliTPCLoader
2656 // and one has to decide where to store the TPC parameters
2659 sprintf(paramName,"75x40_100x60_150x60");
2660 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
2662 AliDebugClass(1,Form("TPC parameters %s found.",paramName));
2664 AliWarningClass("TPC parameters not found. Create new (they may be incorrect)");
2665 paramTPC = new AliTPCParamSR;
2669 // the older version of parameters can be accessed with this code.
2670 // In some cases, we have old parameters saved in the file but
2671 // digits were created with new parameters, it can be distinguish
2672 // by the name of TPC TreeD. The code here is just for the case
2673 // we would need to compare with old data, uncomment it if needed.
2675 // char paramName[50];
2676 // sprintf(paramName,"75x40_100x60");
2677 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
2679 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2681 // sprintf(paramName,"75x40_100x60_150x60");
2682 // paramTPC=(AliTPCParam*)in->Get(paramName);
2684 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2686 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
2688 // paramTPC = new AliTPCParamSR;