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"
75 #include "AliTPCDigitizer.h"
76 #include "AliTPCBuffer.h"
77 #include "AliTPCDDLRawData.h"
79 #include "AliRawReader.h"
80 #include "AliTPCRawStream.h"
83 //_____________________________________________________________________________
84 AliTPC::AliTPC():AliDetector(),
101 // Default constructor
105 // fTrackHitsOld = 0;
106 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
107 fHitType = 4; // ROOT containers
109 fHitType = 2; //default CONTAINERS - based on ROOT structure
115 //_____________________________________________________________________________
116 AliTPC::AliTPC(const char *name, const char *title)
117 : AliDetector(name,title),
133 // Standard constructor
137 // Initialise arrays of hits and digits
138 fHits = new TClonesArray("AliTPChit", 176);
139 gAlice->GetMCApp()->AddHitList(fHits);
141 fTrackHits = new AliTPCTrackHitsV2;
142 fTrackHits->SetHitPrecision(0.002);
143 fTrackHits->SetStepPrecision(0.003);
144 fTrackHits->SetMaxDistance(100);
146 //fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
147 //fTrackHitsOld->SetHitPrecision(0.002);
148 //fTrackHitsOld->SetStepPrecision(0.003);
149 //fTrackHitsOld->SetMaxDistance(100);
152 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
153 fHitType = 4; // ROOT containers
163 // Initialise color attributes
164 //PH SetMarkerColor(kYellow);
167 // Set TPC parameters
171 if (!strcmp(title,"Default")) {
172 fTPCParam = new AliTPCParamSR;
174 AliWarning("In Config.C you must set non-default parameters.");
181 //_____________________________________________________________________________
192 delete fTrackHits; //MI 15.09.2000
193 // delete fTrackHitsOld; //MI 10.12.2001
196 delete [] fNoiseTable;
197 delete [] fActiveSectors;
201 //_____________________________________________________________________________
202 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
205 // Add a hit to the list
208 TClonesArray &lhits = *fHits;
209 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
212 AddHit2(track,vol,hits);
215 //_____________________________________________________________________________
216 void AliTPC::BuildGeometry()
220 // Build TPC ROOT TNode geometry for the event display
225 const int kColorTPC=19;
226 char name[5], title[25];
227 const Double_t kDegrad=TMath::Pi()/180;
228 const Double_t kRaddeg=180./TMath::Pi();
231 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
232 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
234 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
235 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
237 Int_t nLo = fTPCParam->GetNInnerSector()/2;
238 Int_t nHi = fTPCParam->GetNOuterSector()/2;
240 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
241 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
242 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
243 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
246 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
247 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
253 // Get ALICE top node
256 nTop=gAlice->GetGeometry()->GetNode("alice");
260 rl = fTPCParam->GetInnerRadiusLow();
261 ru = fTPCParam->GetInnerRadiusUp();
265 sprintf(name,"LS%2.2d",i);
267 sprintf(title,"TPC low sector %3d",i);
270 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
271 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
272 tubs->SetNumberOfDivisions(1);
274 nNode = new TNode(name,title,name,0,0,0,"");
275 nNode->SetLineColor(kColorTPC);
281 rl = fTPCParam->GetOuterRadiusLow();
282 ru = fTPCParam->GetOuterRadiusUp();
285 sprintf(name,"US%2.2d",i);
287 sprintf(title,"TPC upper sector %d",i);
289 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
290 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
291 tubs->SetNumberOfDivisions(1);
293 nNode = new TNode(name,title,name,0,0,0,"");
294 nNode->SetLineColor(kColorTPC);
300 //_____________________________________________________________________________
301 void AliTPC::CreateMaterials()
303 //-----------------------------------------------
304 // Create Materials for for TPC simulations
305 //-----------------------------------------------
307 //-----------------------------------------------------------------
308 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
309 //-----------------------------------------------------------------
311 Int_t iSXFLD=gAlice->Field()->Integ();
312 Float_t sXMGMX=gAlice->Field()->Max();
314 Float_t amat[5]; // atomic numbers
315 Float_t zmat[5]; // z
316 Float_t wmat[5]; // proportions
322 //***************** Gases *************************
325 //--------------------------------------------------------------
326 // gases - air and CO2
327 //--------------------------------------------------------------
343 AliMixture(10,"CO2",amat,zmat,density,2,wmat);
358 AliMixture(11,"Air",amat,zmat,density,2,wmat);
360 //----------------------------------------------------------------
362 //----------------------------------------------------------------
365 // Drift gases 1 - nonsensitive, 2 - sensitive
366 // Ne-CO2-N (85-10-5)
385 AliMixture(12,"Ne-CO2-N-1",amat,zmat,density,4,wmat);
386 AliMixture(13,"Ne-CO2-N-2",amat,zmat,density,4,wmat);
388 //----------------------------------------------------------------------
390 //----------------------------------------------------------------------
412 AliMixture(14,"Kevlar",amat,zmat,density,-4,wmat);
433 AliMixture(15,"NOMEX",amat,zmat,density,-4,wmat);
451 AliMixture(16,"Makrolon",amat,zmat,density,-3,wmat);
469 AliMixture(17, "Tedlar",amat,zmat,density,-3,wmat);
487 AliMixture(18, "Mylar",amat,zmat,density,-3,wmat);
488 // material for "prepregs"
489 // Epoxy - C14 H20 O3
492 // prepreg1 60% C-fiber, 40% epoxy (vol)
507 AliMixture(19, "Prepreg1",amat,zmat,density,3,wmat);
509 //prepreg2 60% glass-fiber, 40% epoxy
528 AliMixture(20, "Prepreg2",amat,zmat,density,4,wmat);
530 //prepreg3 50% glass-fiber, 50% epoxy
549 AliMixture(21, "Prepreg3",amat,zmat,density,4,wmat);
551 // G10 60% SiO2 40% epoxy
570 AliMixture(22, "G10",amat,zmat,density,4,wmat);
579 AliMaterial(23,"Al",amat[0],zmat[0],density,999.,999.);
581 // Si (for electronics
588 AliMaterial(24,"Si",amat[0],zmat[0],density,999.,999.);
597 AliMaterial(25,"Cu",amat[0],zmat[0],density,999.,999.);
599 // Epoxy - C14 H20 O3
615 AliMixture(26,"Epoxy",amat,zmat,density,-3,wmat);
633 AliMixture(27,"Plexiglas",amat,zmat,density,-3,wmat);
641 AliMaterial(28,"C",amat[0],zmat[0],density,999.,999.);
643 // Fe (steel for the inner heat screen)
651 AliMaterial(29,"Fe",amat[0],zmat[0],density,999.,999.);
653 //----------------------------------------------------------
654 // tracking media for gases
655 //----------------------------------------------------------
657 AliMedium(0, "Air", 11, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
658 AliMedium(1, "Ne-CO2-N-1", 12, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
659 AliMedium(2, "Ne-CO2-N-2", 13, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
660 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
661 AliMedium(20, "Ne-CO2-N-3", 13, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
662 //-----------------------------------------------------------
663 // tracking media for solids
664 //-----------------------------------------------------------
666 AliMedium(4,"Al",23,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
667 AliMedium(5,"Kevlar",14,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
668 AliMedium(6,"Nomex",15,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
669 AliMedium(7,"Makrolon",16,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
670 AliMedium(8,"Mylar",18,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
671 AliMedium(9,"Tedlar",17,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
673 AliMedium(10,"Prepreg1",19,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
674 AliMedium(11,"Prepreg2",20,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
675 AliMedium(12,"Prepreg3",21,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
676 AliMedium(13,"Epoxy",26,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
678 AliMedium(14,"Cu",25,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
679 AliMedium(15,"Si",24,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
680 AliMedium(16,"G10",22,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
681 AliMedium(17,"Plexiglas",27,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
682 AliMedium(18,"Steel",29,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
686 void AliTPC::GenerNoise(Int_t tablesize)
689 //Generate table with noise
696 if (fNoiseTable) delete[] fNoiseTable;
697 fNoiseTable = new Float_t[tablesize];
698 fNoiseDepth = tablesize;
699 fCurrentNoise =0; //!index of the noise in the noise table
701 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
702 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
705 Float_t AliTPC::GetNoise()
707 // get noise from table
708 // if ((fCurrentNoise%10)==0)
709 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
710 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
711 return fNoiseTable[fCurrentNoise++];
712 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
716 Bool_t AliTPC::IsSectorActive(Int_t sec) const
719 // check if the sector is active
720 if (!fActiveSectors) return kTRUE;
721 else return fActiveSectors[sec];
724 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
726 // activate interesting sectors
727 SetTreeAddress();//just for security
728 if (fActiveSectors) delete [] fActiveSectors;
729 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
730 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
731 for (Int_t i=0;i<n;i++)
732 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
736 void AliTPC::SetActiveSectors(Int_t flag)
739 // activate sectors which were hitted by tracks
741 SetTreeAddress();//just for security
742 if (fHitType==0) return; // if Clones hit - not short volume ID information
743 if (fActiveSectors) delete [] fActiveSectors;
744 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
746 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
749 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
753 AliFatal("Can not find TreeH in folder");
756 if (fHitType>1) branch = TreeH()->GetBranch("TPC2");
757 else branch = TreeH()->GetBranch("TPC");
758 Stat_t ntracks = TreeH()->GetEntries();
759 // loop over all hits
760 AliDebug(1,Form("Got %d tracks",ntracks));
762 for(Int_t track=0;track<ntracks;track++) {
765 if (fTrackHits && fHitType&4) {
766 TBranch * br1 = TreeH()->GetBranch("fVolumes");
767 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
768 br1->GetEvent(track);
769 br2->GetEvent(track);
770 Int_t *volumes = fTrackHits->GetVolumes();
771 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
772 fActiveSectors[volumes[j]]=kTRUE;
776 // if (fTrackHitsOld && fHitType&2) {
777 // TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
778 // br->GetEvent(track);
779 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
780 // for (UInt_t j=0;j<ar->GetSize();j++){
781 // fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
790 //_____________________________________________________________________________
791 void AliTPC::Digits2Raw()
793 // convert digits of the current event to raw data
795 static const Int_t kThreshold = 0;
797 fLoader->LoadDigits();
798 TTree* digits = fLoader->TreeD();
800 AliError("No digits tree");
805 AliSimDigits* digrow = &digarr;
806 digits->GetBranch("Segment")->SetAddress(&digrow);
808 const char* fileName = "AliTPCDDL.dat";
809 AliTPCBuffer* buffer = new AliTPCBuffer(fileName);
813 // 2: txt files with digits
814 //BE CAREFUL, verbose level 2 MUST be used only for debugging and
815 //it is highly suggested to use this mode only for debugging digits files
816 //reasonably small, because otherwise the size of the txt files can reach
817 //quickly several MB wasting time and disk space.
818 buffer->SetVerbose(0);
820 Int_t nEntries = Int_t(digits->GetEntries());
821 Int_t previousSector = -1;
823 for (Int_t i = 0; i < nEntries; i++) {
826 fTPCParam->AdjustSectorRow(digarr.GetID(), sector, row);
827 if(previousSector != sector) {
829 previousSector = sector;
832 if (sector < 36) { //inner sector [0;35]
834 //the whole row is written into the output file
835 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
836 sector, subSector, row);
838 //only the pads in the range [37;48] are written into the output file
839 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 1,
840 sector, subSector, row);
842 //only the pads outside the range [37;48] are written into the output file
843 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 2,
844 sector, subSector, row);
847 } else { //outer sector [36;71]
848 if (row == 54) subSector = 2;
849 if ((row != 27) && (row != 76)) {
850 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
851 sector, subSector, row);
852 } else if (row == 27) {
853 //only the pads outside the range [43;46] are written into the output file
854 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 2,
855 sector, subSector, row);
857 //only the pads in the range [43;46] are written into the output file
858 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 1,
859 sector, subSector, row);
860 } else if (row == 76) {
861 //only the pads outside the range [33;88] are written into the output file
862 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 2,
863 sector, subSector, row);
865 //only the pads in the range [33;88] are written into the output file
866 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 1,
867 sector, subSector, row);
873 fLoader->UnloadDigits();
875 AliTPCDDLRawData rawWriter;
876 rawWriter.SetVerbose(0);
878 rawWriter.RawData(fileName);
879 gSystem->Unlink(fileName);
884 //_____________________________________________________________________________
885 Bool_t AliTPC::Raw2SDigits(AliRawReader* rawReader){
886 // Converts the TPC raw data into summable digits
887 // The method is used for merging simulated and
889 if (fLoader->TreeS() == 0x0 ) {
890 fLoader->MakeTree("S");
893 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
895 //setup TPCDigitsArray
896 if(GetDigitsArray()) delete GetDigitsArray();
898 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
899 arr->SetClass("AliSimDigits");
900 arr->Setup(fTPCParam);
901 arr->MakeTree(fLoader->TreeS());
905 // set zero suppression to "0"
906 fTPCParam->SetZeroSup(0);
909 const Int_t kmaxTime = fTPCParam->GetMaxTBin();
910 const Int_t kNIS = fTPCParam->GetNInnerSector();
911 const Int_t kNOS = fTPCParam->GetNOuterSector();
912 const Int_t kNS = kNIS + kNOS;
914 Short_t** allBins = NULL; //array which contains the data for one sector
916 for(Int_t iSector = 0; iSector < kNS; iSector++) {
918 Int_t nRows = fTPCParam->GetNRow(iSector);
919 Int_t nDDLs = 0, indexDDL = 0;
920 if (iSector < kNIS) {
922 indexDDL = iSector * 2;
926 indexDDL = (iSector-kNIS) * 4 + kNIS * 2;
929 // Loas the raw data for corresponding DDLs
931 AliTPCRawStream input(rawReader);
932 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
934 // Alocate and init the array with the sector data
935 allBins = new Short_t*[nRows];
936 for (Int_t iRow = 0; iRow < nRows; iRow++) {
937 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
938 Int_t maxBin = kmaxTime*maxPad;
939 allBins[iRow] = new Short_t[maxBin];
940 memset(allBins[iRow],0,sizeof(Short_t)*maxBin);
943 // Begin loop over altro data
944 while (input.Next()) {
946 if (input.GetSector() != iSector)
947 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",iSector,input.GetSector()));
949 Int_t iRow = input.GetRow();
950 if (iRow < 0 || iRow >= nRows)
951 AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
953 Int_t iPad = input.GetPad();
955 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
957 if (iPad < 0 || iPad >= maxPad)
958 AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
959 iPad, 0, maxPad -1));
961 Int_t iTimeBin = input.GetTime();
962 if ( iTimeBin < 0 || iTimeBin >= kmaxTime)
963 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
964 iTimeBin, 0, kmaxTime -1));
966 Int_t maxBin = kmaxTime*maxPad;
968 if (((iPad*kmaxTime+iTimeBin) >= maxBin) ||
969 ((iPad*kmaxTime+iTimeBin) < 0))
970 AliFatal(Form("Index outside the allowed range"
971 " Sector=%d Row=%d Pad=%d Timebin=%d"
972 " (Max.index=%d)",iSector,iRow,iPad,iTimeBin,maxBin));
974 allBins[iRow][iPad*kmaxTime+iTimeBin] = input.GetSignal();
976 } // End loop over altro data
978 // Now fill the digits array
979 if (fDigitsArray->GetTree()==0) {
980 AliFatal("Tree not set in fDigitsArray");
983 for (Int_t iRow = 0; iRow < nRows; iRow++) {
984 AliDigits * dig = fDigitsArray->CreateRow(iSector,iRow);
986 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
987 for(Int_t iPad = 0; iPad < maxPad; iPad++) {
988 for(Int_t iTimeBin = 0; iTimeBin < kmaxTime; iTimeBin++) {
989 Short_t q = allBins[iRow][iPad*kmaxTime + iTimeBin];
990 if (q <= 0) continue;
992 dig->SetDigitFast((Short_t)q,iTimeBin,iPad);
995 fDigitsArray->StoreRow(iSector,iRow);
996 Int_t ndig = dig->GetDigitSize();
999 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1000 iSector,iRow,ndig));
1002 fDigitsArray->ClearRow(iSector,iRow);
1004 } // end of the sector digitization
1006 for (Int_t iRow = 0; iRow < nRows; iRow++)
1007 delete [] allBins[iRow];
1013 fLoader->WriteSDigits("OVERWRITE");
1015 if(GetDigitsArray()) delete GetDigitsArray();
1016 SetDigitsArray(0x0);
1021 //______________________________________________________________________
1022 AliDigitizer* AliTPC::CreateDigitizer(AliRunDigitizer* manager) const
1024 return new AliTPCDigitizer(manager);
1027 void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
1029 //create digits from summable digits
1030 GenerNoise(500000); //create teble with noise
1032 //conect tree with sSDigits
1033 TTree *t = fLoader->TreeS();
1036 fLoader->LoadSDigits("READ");
1037 t = fLoader->TreeS();
1039 AliError("Can not get input TreeS");
1044 if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
1046 AliSimDigits digarr, *dummy=&digarr;
1047 TBranch* sdb = t->GetBranch("Segment");
1049 AliError("Can not find branch with segments in TreeS.");
1053 sdb->SetAddress(&dummy);
1055 Stat_t nentries = t->GetEntries();
1057 // set zero suppression
1059 fTPCParam->SetZeroSup(2);
1061 // get zero suppression
1063 Int_t zerosup = fTPCParam->GetZeroSup();
1065 //make tree with digits
1067 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1068 arr->SetClass("AliSimDigits");
1069 arr->Setup(fTPCParam);
1070 arr->MakeTree(fLoader->TreeD());
1072 AliTPCParam * par = fTPCParam;
1074 //Loop over segments of the TPC
1076 for (Int_t n=0; n<nentries; n++) {
1079 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1080 AliWarning(Form("Invalid segment ID ! %d",digarr.GetID()));
1083 if (!IsSectorActive(sec)) continue;
1085 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1086 Int_t nrows = digrow->GetNRows();
1087 Int_t ncols = digrow->GetNCols();
1089 digrow->ExpandBuffer();
1090 digarr.ExpandBuffer();
1091 digrow->ExpandTrackBuffer();
1092 digarr.ExpandTrackBuffer();
1095 Short_t * pamp0 = digarr.GetDigits();
1096 Int_t * ptracks0 = digarr.GetTracks();
1097 Short_t * pamp1 = digrow->GetDigits();
1098 Int_t * ptracks1 = digrow->GetTracks();
1099 Int_t nelems =nrows*ncols;
1100 Int_t saturation = fTPCParam->GetADCSat() - 1;
1101 //use internal structure of the AliDigits - for speed reason
1102 //if you cahnge implementation
1103 //of the Alidigits - it must be rewriten -
1104 for (Int_t i= 0; i<nelems; i++){
1105 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1107 if (q>saturation) q=saturation;
1110 ptracks1[0]=ptracks0[0];
1111 ptracks1[nelems]=ptracks0[nelems];
1112 ptracks1[2*nelems]=ptracks0[2*nelems];
1120 arr->StoreRow(sec,row);
1121 arr->ClearRow(sec,row);
1126 fLoader->WriteDigits("OVERWRITE");
1130 //__________________________________________________________________
1131 void AliTPC::SetDefaults(){
1133 // setting the defaults
1136 // Set response functions
1139 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1141 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1143 AliInfo("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...");
1145 param = new AliTPCParamSR();
1148 param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1151 AliFatal("No TPC parameters found");
1155 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1156 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1157 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1158 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1159 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1160 rf->SetOffset(3*param->GetZSigma());
1163 TDirectory *savedir=gDirectory;
1164 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1166 AliFatal("Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !");
1169 prfinner->Read("prf_07504_Gati_056068_d02");
1170 //PH Set different names
1171 s=prfinner->GetGRF()->GetName();
1173 prfinner->GetGRF()->SetName(s.Data());
1175 prfouter1->Read("prf_10006_Gati_047051_d03");
1176 s=prfouter1->GetGRF()->GetName();
1178 prfouter1->GetGRF()->SetName(s.Data());
1180 prfouter2->Read("prf_15006_Gati_047051_d03");
1181 s=prfouter2->GetGRF()->GetName();
1183 prfouter2->GetGRF()->SetName(s.Data());
1188 param->SetInnerPRF(prfinner);
1189 param->SetOuter1PRF(prfouter1);
1190 param->SetOuter2PRF(prfouter2);
1191 param->SetTimeRF(rf);
1201 //__________________________________________________________________
1202 void AliTPC::Hits2Digits()
1205 // creates digits from hits
1207 if (!fTPCParam->IsGeoRead()){
1209 // read transformation matrices for gGeoManager
1211 fTPCParam->ReadGeoMatrices();
1214 fLoader->LoadHits("read");
1215 fLoader->LoadDigits("recreate");
1216 AliRunLoader* runLoader = fLoader->GetRunLoader();
1218 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1219 runLoader->GetEvent(iEvent);
1221 Hits2Digits(iEvent);
1224 fLoader->UnloadHits();
1225 fLoader->UnloadDigits();
1227 //__________________________________________________________________
1228 void AliTPC::Hits2Digits(Int_t eventnumber)
1230 //----------------------------------------------------
1231 // Loop over all sectors for a single event
1232 //----------------------------------------------------
1233 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1234 rl->GetEvent(eventnumber);
1235 if (fLoader->TreeH() == 0x0) {
1236 if(fLoader->LoadHits()) {
1237 AliError("Can not load hits.");
1242 if (fLoader->TreeD() == 0x0 ) {
1243 fLoader->MakeTree("D");
1244 if (fLoader->TreeD() == 0x0 ) {
1245 AliError("Can not get TreeD");
1250 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1251 GenerNoise(500000); //create teble with noise
1253 //setup TPCDigitsArray
1255 if(GetDigitsArray()) delete GetDigitsArray();
1257 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1258 arr->SetClass("AliSimDigits");
1259 arr->Setup(fTPCParam);
1261 arr->MakeTree(fLoader->TreeD());
1262 SetDigitsArray(arr);
1264 fDigitsSwitch=0; // standard digits
1266 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1267 if (IsSectorActive(isec)) {
1268 AliDebug(1,Form("Hits2Digits","Sector %d is active.",isec));
1269 Hits2DigitsSector(isec);
1272 AliDebug(1,Form("Hits2Digits","Sector %d is NOT active.",isec));
1275 fLoader->WriteDigits("OVERWRITE");
1277 //this line prevents the crash in the similar one
1278 //on the beginning of this method
1279 //destructor attempts to reset the tree, which is deleted by the loader
1280 //need to be redesign
1281 if(GetDigitsArray()) delete GetDigitsArray();
1282 SetDigitsArray(0x0);
1286 //__________________________________________________________________
1287 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1290 //-----------------------------------------------------------
1291 // summable digits - 16 bit "ADC", no noise, no saturation
1292 //-----------------------------------------------------------
1294 //----------------------------------------------------
1295 // Loop over all sectors for a single event
1296 //----------------------------------------------------
1298 AliRunLoader* rl = fLoader->GetRunLoader();
1300 rl->GetEvent(eventnumber);
1301 if (fLoader->TreeH() == 0x0) {
1302 if(fLoader->LoadHits()) {
1303 AliError("Can not load hits.");
1310 if (fLoader->TreeS() == 0x0 ) {
1311 fLoader->MakeTree("S");
1314 if(fDefaults == 0) SetDefaults();
1316 GenerNoise(500000); //create table with noise
1317 //setup TPCDigitsArray
1319 if(GetDigitsArray()) delete GetDigitsArray();
1322 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1323 arr->SetClass("AliSimDigits");
1324 arr->Setup(fTPCParam);
1325 arr->MakeTree(fLoader->TreeS());
1327 SetDigitsArray(arr);
1329 fDigitsSwitch=1; // summable digits
1331 // set zero suppression to "0"
1333 fTPCParam->SetZeroSup(0);
1335 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1336 if (IsSectorActive(isec)) {
1337 Hits2DigitsSector(isec);
1340 fLoader->WriteSDigits("OVERWRITE");
1342 //this line prevents the crash in the similar one
1343 //on the beginning of this method
1344 //destructor attempts to reset the tree, which is deleted by the loader
1345 //need to be redesign
1346 if(GetDigitsArray()) delete GetDigitsArray();
1347 SetDigitsArray(0x0);
1349 //__________________________________________________________________
1351 void AliTPC::Hits2SDigits()
1354 //-----------------------------------------------------------
1355 // summable digits - 16 bit "ADC", no noise, no saturation
1356 //-----------------------------------------------------------
1358 if (!fTPCParam->IsGeoRead()){
1360 // read transformation matrices for gGeoManager
1362 fTPCParam->ReadGeoMatrices();
1365 fLoader->LoadHits("read");
1366 fLoader->LoadSDigits("recreate");
1367 AliRunLoader* runLoader = fLoader->GetRunLoader();
1369 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1370 runLoader->GetEvent(iEvent);
1373 Hits2SDigits2(iEvent);
1376 fLoader->UnloadHits();
1377 fLoader->UnloadSDigits();
1379 //_____________________________________________________________________________
1381 void AliTPC::Hits2DigitsSector(Int_t isec)
1383 //-------------------------------------------------------------------
1384 // TPC conversion from hits to digits.
1385 //-------------------------------------------------------------------
1387 //-----------------------------------------------------------------
1388 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1389 //-----------------------------------------------------------------
1391 //-------------------------------------------------------
1392 // Get the access to the track hits
1393 //-------------------------------------------------------
1395 // check if the parameters are set - important if one calls this method
1396 // directly, not from the Hits2Digits
1398 if(fDefaults == 0) SetDefaults();
1400 TTree *tH = TreeH(); // pointer to the hits tree
1402 AliFatal("Can not find TreeH in folder");
1406 Stat_t ntracks = tH->GetEntries();
1410 //-------------------------------------------
1411 // Only if there are any tracks...
1412 //-------------------------------------------
1416 Int_t nrows =fTPCParam->GetNRow(isec);
1418 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1420 MakeSector(isec,nrows,tH,ntracks,row);
1422 //--------------------------------------------------------
1423 // Digitize this sector, row by row
1424 // row[i] is the pointer to the TObjArray of TVectors,
1425 // each one containing electrons accepted on this
1426 // row, assigned into tracks
1427 //--------------------------------------------------------
1431 if (fDigitsArray->GetTree()==0) {
1432 AliFatal("Tree not set in fDigitsArray");
1435 for (i=0;i<nrows;i++){
1437 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1439 DigitizeRow(i,isec,row);
1441 fDigitsArray->StoreRow(isec,i);
1443 Int_t ndig = dig->GetDigitSize();
1446 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1449 fDigitsArray->ClearRow(isec,i);
1452 } // end of the sector digitization
1454 for(i=0;i<nrows+2;i++){
1459 delete [] row; // delete the array of pointers to TObjArray-s
1463 } // end of Hits2DigitsSector
1466 //_____________________________________________________________________________
1467 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1469 //-----------------------------------------------------------
1470 // Single row digitization, coupling from the neighbouring
1471 // rows taken into account
1472 //-----------------------------------------------------------
1474 //-----------------------------------------------------------------
1475 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1476 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1477 //-----------------------------------------------------------------
1479 Float_t zerosup = fTPCParam->GetZeroSup();
1481 fCurrentIndex[1]= isec;
1484 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1485 Int_t nofTbins = fTPCParam->GetMaxTBin();
1486 Int_t indexRange[4];
1488 // Integrated signal for this row
1489 // and a single track signal
1492 TMatrixF *m1 = new TMatrixF(0,nofPads,0,nofTbins); // integrated
1493 TMatrixF *m2 = new TMatrixF(0,nofPads,0,nofTbins); // single
1495 TMatrixF &total = *m1;
1497 // Array of pointers to the label-signal list
1499 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1500 Float_t **pList = new Float_t* [nofDigits];
1504 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1510 for (Int_t row= row1;row<=row2;row++){
1511 Int_t nTracks= rows[row]->GetEntries();
1512 for (i1=0;i1<nTracks;i1++){
1513 fCurrentIndex[2]= row;
1514 fCurrentIndex[3]=irow+1;
1516 m2->Zero(); // clear single track signal matrix
1517 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1518 GetList(trackLabel,nofPads,m2,indexRange,pList);
1520 else GetSignal(rows[row],i1,0,m1,indexRange);
1526 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1528 Float_t fzerosup = zerosup+0.5;
1529 for(Int_t it=0;it<nofTbins;it++){
1530 for(Int_t ip=0;ip<nofPads;ip++){
1532 Float_t q=total(ip,it);
1533 if(fDigitsSwitch == 0){
1535 if(q <=fzerosup) continue; // do not fill zeros
1537 if(q >= fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat() - 1; // saturation
1542 if(q <= 0.) continue; // do not fill zeros
1543 if(q>2000.) q=2000.;
1549 // "real" signal or electronic noise (list = -1)?
1552 for(Int_t j1=0;j1<3;j1++){
1553 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1558 <A NAME="AliDigits"></A>
1559 using of AliDigits object
1562 dig->SetDigitFast((Short_t)q,it,ip);
1563 if (fDigitsArray->IsSimulated()) {
1564 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1565 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1566 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1569 } // end of loop over time buckets
1570 } // end of lop over pads
1573 // This row has been digitized, delete nonused stuff
1576 for(lp=0;lp<nofDigits;lp++){
1577 if(pList[lp]) delete [] pList[lp];
1585 } // end of DigitizeRow
1587 //_____________________________________________________________________________
1589 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1590 TMatrixF *m1, TMatrixF *m2,Int_t *indexRange)
1593 //---------------------------------------------------------------
1594 // Calculates 2-D signal (pad,time) for a single track,
1595 // returns a pointer to the signal matrix and the track label
1596 // No digitization is performed at this level!!!
1597 //---------------------------------------------------------------
1599 //-----------------------------------------------------------------
1600 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1601 // Modified: Marian Ivanov
1602 //-----------------------------------------------------------------
1606 tv = (TVector*)p1->At(ntr); // pointer to a track
1609 Float_t label = v(0);
1610 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1))/2;
1612 Int_t nElectrons = (tv->GetNrows()-1)/5;
1613 indexRange[0]=9999; // min pad
1614 indexRange[1]=-1; // max pad
1615 indexRange[2]=9999; //min time
1616 indexRange[3]=-1; // max time
1618 TMatrixF &signal = *m1;
1619 TMatrixF &total = *m2;
1621 // Loop over all electrons
1623 for(Int_t nel=0; nel<nElectrons; nel++){
1625 Float_t aval = v(idx+4);
1626 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1627 Float_t xyz[4]={v(idx+1),v(idx+2),v(idx+3),v(idx+5)};
1628 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
1630 Int_t *index = fTPCParam->GetResBin(0);
1631 Float_t *weight = & (fTPCParam->GetResWeight(0));
1633 if (n>0) for (Int_t i =0; i<n; i++){
1634 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1637 Int_t time=index[2];
1638 Float_t qweight = *(weight)*eltoadcfac;
1640 if (m1!=0) signal(pad,time)+=qweight;
1641 total(pad,time)+=qweight;
1642 if (indexRange[0]>pad) indexRange[0]=pad;
1643 if (indexRange[1]<pad) indexRange[1]=pad;
1644 if (indexRange[2]>time) indexRange[2]=time;
1645 if (indexRange[3]<time) indexRange[3]=time;
1652 } // end of loop over electrons
1654 return label; // returns track label when finished
1657 //_____________________________________________________________________________
1658 void AliTPC::GetList(Float_t label,Int_t np,TMatrixF *m,
1659 Int_t *indexRange, Float_t **pList)
1661 //----------------------------------------------------------------------
1662 // Updates the list of tracks contributing to digits for a given row
1663 //----------------------------------------------------------------------
1665 //-----------------------------------------------------------------
1666 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1667 //-----------------------------------------------------------------
1669 TMatrixF &signal = *m;
1671 // lop over nonzero digits
1673 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1674 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1677 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1679 if(signal(ip,it)<0.5) continue;
1681 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1683 if(!pList[globalIndex]){
1686 // Create new list (6 elements - 3 signals and 3 labels),
1689 pList[globalIndex] = new Float_t [6];
1693 *pList[globalIndex] = -1.;
1694 *(pList[globalIndex]+1) = -1.;
1695 *(pList[globalIndex]+2) = -1.;
1696 *(pList[globalIndex]+3) = -1.;
1697 *(pList[globalIndex]+4) = -1.;
1698 *(pList[globalIndex]+5) = -1.;
1700 *pList[globalIndex] = label;
1701 *(pList[globalIndex]+3) = signal(ip,it);
1705 // check the signal magnitude
1707 Float_t highest = *(pList[globalIndex]+3);
1708 Float_t middle = *(pList[globalIndex]+4);
1709 Float_t lowest = *(pList[globalIndex]+5);
1712 // compare the new signal with already existing list
1715 if(signal(ip,it)<lowest) continue; // neglect this track
1719 if (signal(ip,it)>highest){
1720 *(pList[globalIndex]+5) = middle;
1721 *(pList[globalIndex]+4) = highest;
1722 *(pList[globalIndex]+3) = signal(ip,it);
1724 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1725 *(pList[globalIndex]+1) = *pList[globalIndex];
1726 *pList[globalIndex] = label;
1728 else if (signal(ip,it)>middle){
1729 *(pList[globalIndex]+5) = middle;
1730 *(pList[globalIndex]+4) = signal(ip,it);
1732 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1733 *(pList[globalIndex]+1) = label;
1736 *(pList[globalIndex]+5) = signal(ip,it);
1737 *(pList[globalIndex]+2) = label;
1741 } // end of loop over pads
1742 } // end of loop over time bins
1745 //___________________________________________________________________
1746 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1747 Stat_t ntracks,TObjArray **row)
1750 //-----------------------------------------------------------------
1751 // Prepares the sector digitization, creates the vectors of
1752 // tracks for each row of this sector. The track vector
1753 // contains the track label and the position of electrons.
1754 //-----------------------------------------------------------------
1756 //-----------------------------------------------------------------
1757 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1758 //-----------------------------------------------------------------
1760 Float_t gasgain = fTPCParam->GetGasGain();
1764 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1767 if (fHitType>1) branch = TH->GetBranch("TPC2");
1768 else branch = TH->GetBranch("TPC");
1771 //----------------------------------------------
1772 // Create TObjArray-s, one for each row,
1773 // each TObjArray will store the TVectors
1774 // of electrons, one TVectors per each track.
1775 //----------------------------------------------
1777 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
1778 TVector **tracks = new TVector* [nrows+2]; //pointers to the track vectors
1780 for(i=0; i<nrows+2; i++){
1781 row[i] = new TObjArray;
1788 //--------------------------------------------------------------------
1789 // Loop over tracks, the "track" contains the full history
1790 //--------------------------------------------------------------------
1792 Int_t previousTrack,currentTrack;
1793 previousTrack = -1; // nothing to store so far!
1795 for(Int_t track=0;track<ntracks;track++){
1796 Bool_t isInSector=kTRUE;
1798 isInSector = TrackInVolume(isec,track);
1799 if (!isInSector) continue;
1801 branch->GetEntry(track); // get next track
1805 tpcHit = (AliTPChit*)FirstHit(-1);
1807 //--------------------------------------------------------------
1809 //--------------------------------------------------------------
1814 Int_t sector=tpcHit->fSector; // sector number
1816 tpcHit = (AliTPChit*) NextHit();
1820 // Remove hits which arrive before the TPC opening gate signal
1821 if(((fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
1822 /fTPCParam->GetDriftV()+tpcHit->Time())<fTPCParam->GetGateDelay()) {
1823 tpcHit = (AliTPChit*) NextHit();
1827 currentTrack = tpcHit->Track(); // track number
1829 if(currentTrack != previousTrack){
1831 // store already filled fTrack
1833 for(i=0;i<nrows+2;i++){
1834 if(previousTrack != -1){
1835 if(nofElectrons[i]>0){
1836 TVector &v = *tracks[i];
1837 v(0) = previousTrack;
1838 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
1839 row[i]->Add(tracks[i]);
1842 delete tracks[i]; // delete empty TVector
1848 tracks[i] = new TVector(601); // TVectors for the next fTrack
1850 } // end of loop over rows
1852 previousTrack=currentTrack; // update track label
1855 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1857 //---------------------------------------------------
1858 // Calculate the electron attachment probability
1859 //---------------------------------------------------
1862 Float_t time = 1.e6*(fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
1863 /fTPCParam->GetDriftV();
1865 Float_t attProb = fTPCParam->GetAttCoef()*
1866 fTPCParam->GetOxyCont()*time; // fraction!
1868 //-----------------------------------------------
1869 // Loop over electrons
1870 //-----------------------------------------------
1873 for(Int_t nel=0;nel<qI;nel++){
1874 // skip if electron lost due to the attachment
1875 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1880 // protection for the nonphysical avalanche size (10**6 maximum)
1882 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
1883 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
1886 TransportElectron(xyz,index);
1888 fTPCParam->GetPadRow(xyz,index);
1889 // Electron track time (for pileup simulation)
1890 xyz[4] = tpcHit->Time()/fTPCParam->GetTSample();
1891 // row 0 - cross talk from the innermost row
1892 // row fNRow+1 cross talk from the outermost row
1893 rowNumber = index[2]+1;
1894 //transform position to local digit coordinates
1895 //relative to nearest pad row
1896 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
1898 if (isec <fTPCParam->GetNInnerSector()) {
1899 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
1900 y1 = fTPCParam->GetYInner(rowNumber);
1903 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
1904 y1 = fTPCParam->GetYOuter(rowNumber);
1906 // gain inefficiency at the wires edges - linear
1909 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));
1911 nofElectrons[rowNumber]++;
1912 //----------------------------------
1913 // Expand vector if necessary
1914 //----------------------------------
1915 if(nofElectrons[rowNumber]>120){
1916 Int_t range = tracks[rowNumber]->GetNrows();
1917 if((nofElectrons[rowNumber])>(range-1)/5){
1919 tracks[rowNumber]->ResizeTo(range+500); // Add 100 electrons
1923 TVector &v = *tracks[rowNumber];
1924 Int_t idx = 5*nofElectrons[rowNumber]-4;
1925 Real_t * position = &(((TVector&)v)(idx)); //make code faster
1926 memcpy(position,xyz,5*sizeof(Float_t));
1928 } // end of loop over electrons
1930 tpcHit = (AliTPChit*)NextHit();
1932 } // end of loop over hits
1933 } // end of loop over tracks
1936 // store remaining track (the last one) if not empty
1939 for(i=0;i<nrows+2;i++){
1940 if(nofElectrons[i]>0){
1941 TVector &v = *tracks[i];
1942 v(0) = previousTrack;
1943 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
1944 row[i]->Add(tracks[i]);
1953 delete [] nofElectrons;
1955 } // end of MakeSector
1958 //_____________________________________________________________________________
1962 // Initialise TPC detector after definition of geometry
1964 AliDebug(1,"*********************************************");
1967 //_____________________________________________________________________________
1968 void AliTPC::MakeBranch(Option_t* option)
1971 // Create Tree branches for the TPC.
1974 Int_t buffersize = 4000;
1975 char branchname[10];
1976 sprintf(branchname,"%s",GetName());
1978 const char *h = strstr(option,"H");
1980 if ( h && (fHitType<=1) && (fHits == 0x0)) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
1982 AliDetector::MakeBranch(option);
1984 const char *d = strstr(option,"D");
1986 if (fDigits && fLoader->TreeD() && d) {
1987 MakeBranchInTree(gAlice->TreeD(), branchname, &fDigits, buffersize, 0);
1990 if (fHitType>1) MakeBranch2(option,0); // MI change 14.09.2000
1993 //_____________________________________________________________________________
1994 void AliTPC::ResetDigits()
1997 // Reset number of digits and the digits array for this detector
2000 if (fDigits) fDigits->Clear();
2005 //_____________________________________________________________________________
2006 void AliTPC::SetSens(Int_t sens)
2009 //-------------------------------------------------------------
2010 // Activates/deactivates the sensitive strips at the center of
2011 // the pad row -- this is for the space-point resolution calculations
2012 //-------------------------------------------------------------
2014 //-----------------------------------------------------------------
2015 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2016 //-----------------------------------------------------------------
2022 void AliTPC::SetSide(Float_t side=0.)
2024 // choice of the TPC side
2029 //_____________________________________________________________________________
2031 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2034 // electron transport taking into account:
2036 // 2.ExB at the wires
2037 // 3. nonisochronity
2039 // xyz and index must be already transformed to system 1
2042 fTPCParam->Transform1to2(xyz,index);
2045 Float_t driftl=xyz[2];
2046 if(driftl<0.01) driftl=0.01;
2047 driftl=TMath::Sqrt(driftl);
2048 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2049 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2050 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2051 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2052 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2056 if (fTPCParam->GetMWPCReadout()==kTRUE){
2057 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2058 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2060 //add nonisochronity (not implemented yet)
2064 //______________________________________________________________________
2065 AliTPChit::AliTPChit()
2077 //_____________________________________________________________________________
2078 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits)
2079 :AliHit(shunt,track),
2086 // Creates a TPC hit object
2097 //________________________________________________________________________
2098 // Additional code because of the AliTPCTrackHitsV2
2100 void AliTPC::MakeBranch2(Option_t *option,const char */*file*/)
2103 // Create a new branch in the current Root Tree
2104 // The branch of fHits is automatically split
2105 // MI change 14.09.2000
2107 if (fHitType<2) return;
2108 char branchname[10];
2109 sprintf(branchname,"%s2",GetName());
2111 // Get the pointer to the header
2112 const char *cH = strstr(option,"H");
2114 if (fTrackHits && TreeH() && cH && fHitType&4) {
2115 AliDebug(1,"Making branch for Type 4 Hits");
2116 TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2119 // if (fTrackHitsOld && TreeH() && cH && fHitType&2) {
2120 // AliDebug(1,"Making branch for Type 2 Hits");
2121 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2122 // TreeH(),fBufferSize,99);
2123 // TreeH()->GetListOfBranches()->Add(branch);
2127 void AliTPC::SetTreeAddress()
2129 //Sets tree address for hits
2131 if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2132 AliDetector::SetTreeAddress();
2134 if (fHitType>1) SetTreeAddress2();
2137 void AliTPC::SetTreeAddress2()
2140 // Set branch address for the TrackHits Tree
2145 char branchname[20];
2146 sprintf(branchname,"%s2",GetName());
2148 // Branch address for hit tree
2149 TTree *treeH = TreeH();
2150 if ((treeH)&&(fHitType&4)) {
2151 branch = treeH->GetBranch(branchname);
2153 branch->SetAddress(&fTrackHits);
2154 AliDebug(1,"fHitType&4 Setting");
2157 AliDebug(1,"fHitType&4 Failed (can not find branch)");
2160 // if ((treeH)&&(fHitType&2)) {
2161 // branch = treeH->GetBranch(branchname);
2163 // branch->SetAddress(&fTrackHitsOld);
2164 // AliDebug(1,"fHitType&2 Setting");
2167 // AliDebug(1,"fHitType&2 Failed (can not find branch)");
2169 //set address to TREETR
2171 TTree *treeTR = TreeTR();
2172 if (treeTR && fTrackReferences) {
2173 branch = treeTR->GetBranch(GetName());
2174 if (branch) branch->SetAddress(&fTrackReferences);
2179 void AliTPC::FinishPrimary()
2181 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2182 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2186 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2189 // add hit to the list
2192 int primary = gAlice->GetMCApp()->GetPrimary(track);
2193 gAlice->GetMCApp()->Particle(primary)->SetBit(kKeepBit);
2197 gAlice->GetMCApp()->FlagTrack(track);
2199 if (fTrackHits && fHitType&4)
2200 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2201 hits[1],hits[2],(Int_t)hits[3],hits[4]);
2202 // if (fTrackHitsOld &&fHitType&2 )
2203 // fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2204 // hits[1],hits[2],(Int_t)hits[3]);
2208 void AliTPC::ResetHits()
2210 if (fHitType&1) AliDetector::ResetHits();
2211 if (fHitType>1) ResetHits2();
2214 void AliTPC::ResetHits2()
2218 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2219 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2223 AliHit* AliTPC::FirstHit(Int_t track)
2225 if (fHitType>1) return FirstHit2(track);
2226 return AliDetector::FirstHit(track);
2228 AliHit* AliTPC::NextHit()
2233 if (fHitType>1) return NextHit2();
2235 return AliDetector::NextHit();
2238 AliHit* AliTPC::FirstHit2(Int_t track)
2241 // Initialise the hit iterator
2242 // Return the address of the first hit for track
2243 // If track>=0 the track is read from disk
2244 // while if track<0 the first hit of the current
2245 // track is returned
2248 gAlice->ResetHits();
2249 TreeH()->GetEvent(track);
2252 if (fTrackHits && fHitType&4) {
2253 fTrackHits->First();
2254 return fTrackHits->GetHit();
2256 // if (fTrackHitsOld && fHitType&2) {
2257 // fTrackHitsOld->First();
2258 // return fTrackHitsOld->GetHit();
2264 AliHit* AliTPC::NextHit2()
2267 //Return the next hit for the current track
2270 // if (fTrackHitsOld && fHitType&2) {
2271 // fTrackHitsOld->Next();
2272 // return fTrackHitsOld->GetHit();
2276 return fTrackHits->GetHit();
2282 void AliTPC::LoadPoints(Int_t)
2287 if(fHitType==1) AliDetector::LoadPoints(a);
2288 else LoadPoints2(a);
2292 void AliTPC::RemapTrackHitIDs(Int_t *map)
2297 if (!fTrackHits) return;
2299 // if (fTrackHitsOld && fHitType&2){
2300 // AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2301 // for (UInt_t i=0;i<arr->GetSize();i++){
2302 // AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2303 // info->fTrackID = map[info->fTrackID];
2306 // if (fTrackHitsOld && fHitType&4){
2307 if (fTrackHits && fHitType&4){
2308 TClonesArray * arr = fTrackHits->GetArray();;
2309 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2310 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2311 info->SetTrackID(map[info->GetTrackID()]);
2316 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2318 //return bool information - is track in given volume
2319 //load only part of the track information
2320 //return true if current track is in volume
2323 // if (fTrackHitsOld && fHitType&2) {
2324 // TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
2325 // br->GetEvent(track);
2326 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2327 // for (UInt_t j=0;j<ar->GetSize();j++){
2328 // if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2332 if (fTrackHits && fHitType&4) {
2333 TBranch * br1 = TreeH()->GetBranch("fVolumes");
2334 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
2335 br2->GetEvent(track);
2336 br1->GetEvent(track);
2337 Int_t *volumes = fTrackHits->GetVolumes();
2338 Int_t nvolumes = fTrackHits->GetNVolumes();
2339 if (!volumes && nvolumes>0) {
2340 AliWarning(Form("Problematic track\t%d\t%d",track,nvolumes));
2343 for (Int_t j=0;j<nvolumes; j++)
2344 if (volumes[j]==id) return kTRUE;
2348 TBranch * br = TreeH()->GetBranch("fSector");
2349 br->GetEvent(track);
2350 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2351 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2358 //_____________________________________________________________________________
2359 void AliTPC::LoadPoints2(Int_t)
2362 // Store x, y, z of all hits in memory
2364 // if (fTrackHits == 0 && fTrackHitsOld==0) return;
2365 if (fTrackHits == 0 ) return;
2368 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2369 // if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2371 if (nhits == 0) return;
2372 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2373 if (fPoints == 0) fPoints = new TObjArray(tracks);
2376 Int_t *ntrk=new Int_t[tracks];
2377 Int_t *limi=new Int_t[tracks];
2378 Float_t **coor=new Float_t*[tracks];
2379 for(Int_t i=0;i<tracks;i++) {
2385 AliPoints *points = 0;
2388 Int_t chunk=nhits/4+1;
2390 // Loop over all the hits and store their position
2392 ahit = FirstHit2(-1);
2394 trk=ahit->GetTrack();
2395 if(ntrk[trk]==limi[trk]) {
2397 // Initialise a new track
2398 fp=new Float_t[3*(limi[trk]+chunk)];
2400 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2401 delete [] coor[trk];
2408 fp[3*ntrk[trk] ] = ahit->X();
2409 fp[3*ntrk[trk]+1] = ahit->Y();
2410 fp[3*ntrk[trk]+2] = ahit->Z();
2418 for(trk=0; trk<tracks; ++trk) {
2420 points = new AliPoints();
2421 points->SetMarkerColor(kYellow); //PH kYellow it the default in TPC
2422 points->SetMarkerSize(1);//PH Default size=1
2423 points->SetDetector(this);
2424 points->SetParticle(trk);
2425 points->SetPolyMarker(ntrk[trk],coor[trk],1); // Default style=1
2426 fPoints->AddAt(points,trk);
2427 delete [] coor[trk];
2437 //_____________________________________________________________________________
2438 void AliTPC::LoadPoints3(Int_t)
2441 // Store x, y, z of all hits in memory
2442 // - only intersection point with pad row
2443 if (fTrackHits == 0) return;
2445 Int_t nhits = fTrackHits->GetEntriesFast();
2446 if (nhits == 0) return;
2447 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2448 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2449 fPoints->Expand(2*tracks);
2452 Int_t *ntrk=new Int_t[tracks];
2453 Int_t *limi=new Int_t[tracks];
2454 Float_t **coor=new Float_t*[tracks];
2455 for(Int_t i=0;i<tracks;i++) {
2461 AliPoints *points = 0;
2464 Int_t chunk=nhits/4+1;
2466 // Loop over all the hits and store their position
2468 ahit = FirstHit2(-1);
2472 trk=ahit->GetTrack();
2473 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2474 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2475 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2476 if (currentrow!=lastrow){
2477 lastrow = currentrow;
2478 //later calculate intersection point
2479 if(ntrk[trk]==limi[trk]) {
2481 // Initialise a new track
2482 fp=new Float_t[3*(limi[trk]+chunk)];
2484 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2485 delete [] coor[trk];
2492 fp[3*ntrk[trk] ] = ahit->X();
2493 fp[3*ntrk[trk]+1] = ahit->Y();
2494 fp[3*ntrk[trk]+2] = ahit->Z();
2501 for(trk=0; trk<tracks; ++trk) {
2503 points = new AliPoints();
2504 points->SetMarkerColor(kMagenta); //PH kYellow + 1 is kMagenta
2505 points->SetMarkerStyle(5);
2506 points->SetMarkerSize(0.2);
2507 points->SetDetector(this);
2508 points->SetParticle(trk);
2509 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2510 fPoints->AddAt(points,tracks+trk);
2511 delete [] coor[trk];
2522 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
2525 fLoader = new AliTPCLoader(GetName(),topfoldername);
2529 ////////////////////////////////////////////////////////////////////////
2530 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
2532 // load TPC paarmeters from a given file or create new if the object
2533 // is not found there
2534 // 12/05/2003 This method should be moved to the AliTPCLoader
2535 // and one has to decide where to store the TPC parameters
2538 sprintf(paramName,"75x40_100x60_150x60");
2539 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
2541 AliDebugClass(1,Form("TPC parameters %s found.",paramName));
2543 AliWarningClass("TPC parameters not found. Create new (they may be incorrect)");
2544 paramTPC = new AliTPCParamSR;
2548 // the older version of parameters can be accessed with this code.
2549 // In some cases, we have old parameters saved in the file but
2550 // digits were created with new parameters, it can be distinguish
2551 // by the name of TPC TreeD. The code here is just for the case
2552 // we would need to compare with old data, uncomment it if needed.
2554 // char paramName[50];
2555 // sprintf(paramName,"75x40_100x60");
2556 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
2558 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2560 // sprintf(paramName,"75x40_100x60_150x60");
2561 // paramTPC=(AliTPCParam*)in->Get(paramName);
2563 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2565 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
2567 // paramTPC = new AliTPCParamSR;