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>
41 #include <TGeoGlobalMagField.h>
42 #include <TInterpreter.h>
45 #include <TObjectTable.h>
46 #include <TParticle.h>
49 #include <TStopwatch.h>
54 #include <TVirtualMC.h>
55 #include <TParameter.h>
57 #include "AliDigits.h"
60 #include "AliRunLoader.h"
61 #include "AliSimDigits.h"
64 #include "AliTPCDigitsArray.h"
65 #include "AliTPCLoader.h"
66 #include "AliTPCPRF2D.h"
67 #include "AliTPCParamSR.h"
68 #include "AliTPCRF1D.h"
69 #include "AliTPCTrackHitsV2.h"
70 #include "AliTrackReference.h"
73 #include "AliTPCDigitizer.h"
74 #include "AliTPCBuffer.h"
75 #include "AliTPCDDLRawData.h"
77 #include "AliTPCcalibDB.h"
78 #include "AliTPCCalPad.h"
79 #include "AliTPCCalROC.h"
80 #include "AliTPCExB.h"
81 #include "AliRawReader.h"
82 #include "AliTPCRawStreamV3.h"
83 #include "TTreeStream.h"
86 //_____________________________________________________________________________
87 AliTPC::AliTPC():AliDetector(),
97 fPrimaryIonisation(0),
108 // Default constructor
112 // fTrackHitsOld = 0;
113 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
114 fHitType = 4; // ROOT containers
116 fHitType = 2; //default CONTAINERS - based on ROOT structure
120 //_____________________________________________________________________________
121 AliTPC::AliTPC(const char *name, const char *title)
122 : AliDetector(name,title),
132 fPrimaryIonisation(0),
143 // Standard constructor
147 // Initialise arrays of hits and digits
148 fHits = new TClonesArray("AliTPChit", 176);
149 gAlice->GetMCApp()->AddHitList(fHits);
151 fTrackHits = new AliTPCTrackHitsV2;
152 fTrackHits->SetHitPrecision(0.002);
153 fTrackHits->SetStepPrecision(0.003);
154 fTrackHits->SetMaxDistance(100);
156 //fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
157 //fTrackHitsOld->SetHitPrecision(0.002);
158 //fTrackHitsOld->SetStepPrecision(0.003);
159 //fTrackHitsOld->SetMaxDistance(100);
162 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
163 fHitType = 4; // ROOT containers
173 // Initialise color attributes
174 //PH SetMarkerColor(kYellow);
177 // Set TPC parameters
181 if (!strcmp(title,"Default")) {
182 //fTPCParam = new AliTPCParamSR;
183 fTPCParam = AliTPCcalibDB::Instance()->GetParameters();
185 AliWarning("In Config.C you must set non-default parameters.");
190 void AliTPC::CreateDebugStremer(){
192 // Create Debug streamer to check simulation
194 fDebugStreamer = new TTreeSRedirector("TPCSimdebug.root");
196 //_____________________________________________________________________________
207 delete fTrackHits; //MI 15.09.2000
208 // delete fTrackHitsOld; //MI 10.12.2001
211 delete [] fNoiseTable;
212 delete [] fActiveSectors;
213 if (fDebugStreamer) delete fDebugStreamer;
216 //_____________________________________________________________________________
217 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
220 // Add a hit to the list
223 TClonesArray &lhits = *fHits;
224 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
227 AddHit2(track,vol,hits);
230 //_____________________________________________________________________________
231 void AliTPC::CreateMaterials()
233 //-----------------------------------------------
234 // Create Materials for for TPC simulations
235 //-----------------------------------------------
237 //-----------------------------------------------------------------
238 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
239 //-----------------------------------------------------------------
241 Int_t iSXFLD=((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ();
242 Float_t sXMGMX=((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();
244 Float_t amat[5]; // atomic numbers
245 Float_t zmat[5]; // z
246 Float_t wmat[5]; // proportions
252 //***************** Gases *************************
255 //--------------------------------------------------------------
256 // gases - air and CO2
257 //--------------------------------------------------------------
273 AliMixture(10,"CO2",amat,zmat,density,2,wmat);
288 AliMixture(11,"Air",amat,zmat,density,2,wmat);
290 //----------------------------------------------------------------
292 //----------------------------------------------------------------
295 // Drift gases 1 - nonsensitive, 2 - sensitive
296 // Ne-CO2-N2 (90-10-5) (volume) values at 20deg and 1 atm.
297 // rho(Ne) = 0.839 g/cm^3, rho(CO2) = 1.842 g/cm^3, rho(N2) = 1.165 g/cm^3
298 // for the calculation - everything is normalized to 1
317 AliMixture(12,"Ne-CO2-N-1",amat,zmat,density,4,wmat);
318 AliMixture(13,"Ne-CO2-N-2",amat,zmat,density,4,wmat);
319 AliMixture(35,"Ne-CO2-N-3",amat,zmat,density,4,wmat);
320 //----------------------------------------------------------------------
322 //----------------------------------------------------------------------
344 AliMixture(14,"Kevlar",amat,zmat,density,-4,wmat);
365 AliMixture(15,"NOMEX",amat,zmat,density,-4,wmat);
383 AliMixture(16,"Makrolon",amat,zmat,density,-3,wmat);
401 AliMixture(17, "Tedlar",amat,zmat,density,-3,wmat);
419 AliMixture(18, "Mylar",amat,zmat,density,-3,wmat);
420 // material for "prepregs"
421 // Epoxy - C14 H20 O3
424 // prepreg1 60% C-fiber, 40% epoxy (vol)
439 AliMixture(19, "Prepreg1",amat,zmat,density,3,wmat);
441 //prepreg2 60% glass-fiber, 40% epoxy
460 AliMixture(20, "Prepreg2",amat,zmat,density,4,wmat);
462 //prepreg3 50% glass-fiber, 50% epoxy
481 AliMixture(21, "Prepreg3",amat,zmat,density,4,wmat);
483 // G10 60% SiO2 40% epoxy
502 AliMixture(22, "G10",amat,zmat,density,4,wmat);
511 AliMaterial(23,"Al",amat[0],zmat[0],density,999.,999.);
513 // Si (for electronics
520 AliMaterial(24,"Si",amat[0],zmat[0],density,999.,999.);
529 AliMaterial(25,"Cu",amat[0],zmat[0],density,999.,999.);
547 AliMixture(33,"Brass",amat,zmat,density,2,wmat);
549 // Epoxy - C14 H20 O3
565 AliMixture(26,"Epoxy",amat,zmat,density,-3,wmat);
567 // epoxy film - 90% epoxy, 10% glass fiber
587 AliMixture(34, "Epoxy-film",amat,zmat,density,4,wmat);
605 AliMixture(27,"Plexiglas",amat,zmat,density,-3,wmat);
613 AliMaterial(28,"C",amat[0],zmat[0],density,999.,999.);
615 // Fe (steel for the inner heat screen)
623 AliMaterial(29,"Fe",amat[0],zmat[0],density,999.,999.);
625 // Peek - (C6H4-O-OC6H4-O-C6H4-CO)n
640 AliMixture(30,"Peek",amat,zmat,density,-3,wmat);
655 AliMixture(31,"Alumina",amat,zmat,density,-2,wmat);
674 AliMixture(32,"Water",amat,zmat,density,-2,wmat);
677 //----------------------------------------------------------
678 // tracking media for gases
679 //----------------------------------------------------------
681 AliMedium(0, "Air", 11, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
682 AliMedium(1, "Ne-CO2-N-1", 12, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
683 AliMedium(2, "Ne-CO2-N-2", 13, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
684 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
685 AliMedium(20, "Ne-CO2-N-3", 35, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
686 //-----------------------------------------------------------
687 // tracking media for solids
688 //-----------------------------------------------------------
690 AliMedium(4,"Al",23,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
691 AliMedium(5,"Kevlar",14,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
692 AliMedium(6,"Nomex",15,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
693 AliMedium(7,"Makrolon",16,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
694 AliMedium(8,"Mylar",18,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
695 AliMedium(9,"Tedlar",17,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
697 AliMedium(10,"Prepreg1",19,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
698 AliMedium(11,"Prepreg2",20,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
699 AliMedium(12,"Prepreg3",21,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
700 AliMedium(13,"Epoxy",26,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
702 AliMedium(14,"Cu",25,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
703 AliMedium(15,"Si",24,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
704 AliMedium(16,"G10",22,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
705 AliMedium(17,"Plexiglas",27,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
706 AliMedium(18,"Steel",29,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
707 AliMedium(19,"Peek",30,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
708 AliMedium(21,"Alumina",31,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
709 AliMedium(22,"Water",32,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
710 AliMedium(23,"Brass",33,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
711 AliMedium(24,"Epoxyfm",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
714 void AliTPC::GenerNoise(Int_t tablesize)
717 //Generate table with noise
724 if (fNoiseTable) delete[] fNoiseTable;
725 fNoiseTable = new Float_t[tablesize];
726 fNoiseDepth = tablesize;
727 fCurrentNoise =0; //!index of the noise in the noise table
729 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
730 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
733 Float_t AliTPC::GetNoise()
735 // get noise from table
736 // if ((fCurrentNoise%10)==0)
737 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
738 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
739 return fNoiseTable[fCurrentNoise++];
740 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
744 Bool_t AliTPC::IsSectorActive(Int_t sec) const
747 // check if the sector is active
748 if (!fActiveSectors) return kTRUE;
749 else return fActiveSectors[sec];
752 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
754 // activate interesting sectors
755 SetTreeAddress();//just for security
756 if (!fActiveSectors) fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
757 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
758 for (Int_t i=0;i<n;i++)
759 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
763 void AliTPC::SetActiveSectors(Int_t flag)
766 // activate sectors which were hitted by tracks
768 SetTreeAddress();//just for security
769 if (fHitType==0) return; // if Clones hit - not short volume ID information
770 if (!fActiveSectors) fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
772 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
775 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
777 if (fLoader->TreeH() == 0x0)
779 AliFatal("Can not find TreeH in folder");
782 if (fHitType>1) branch = fLoader->TreeH()->GetBranch("TPC2");
783 else branch = fLoader->TreeH()->GetBranch("TPC");
784 Stat_t ntracks = fLoader->TreeH()->GetEntries();
785 // loop over all hits
786 AliDebug(1,Form("Got %d tracks", (Int_t) ntracks));
788 for(Int_t track=0;track<ntracks;track++) {
791 if (fTrackHits && fHitType&4) {
792 TBranch * br1 = fLoader->TreeH()->GetBranch("fVolumes");
793 TBranch * br2 = fLoader->TreeH()->GetBranch("fNVolumes");
794 br1->GetEvent(track);
795 br2->GetEvent(track);
796 Int_t *volumes = fTrackHits->GetVolumes();
797 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++) {
798 if (volumes[j]>-1 && volumes[j]<fTPCParam->GetNSector()) {
799 fActiveSectors[volumes[j]]=kTRUE;
802 AliError(Form("Volume %d -> sector number %d is outside (0..%d)",
805 fTPCParam->GetNSector()));
811 // if (fTrackHitsOld && fHitType&2) {
812 // TBranch * br = fLoader->TreeH()->GetBranch("fTrackHitsInfo");
813 // br->GetEvent(track);
814 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
815 // for (UInt_t j=0;j<ar->GetSize();j++){
816 // fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
825 //_____________________________________________________________________________
826 void AliTPC::Digits2Raw()
828 // convert digits of the current event to raw data
830 static const Int_t kThreshold = 0;
832 fLoader->LoadDigits();
833 TTree* digits = fLoader->TreeD();
835 AliError("No digits tree");
841 AliSimDigits* digrow = &digarr;
842 digits->GetBranch("Segment")->SetAddress(&digrow);
844 const char* fileName = "AliTPCDDL.dat";
845 AliTPCBuffer* buffer = new AliTPCBuffer(fileName);
849 // 2: txt files with digits
850 //BE CAREFUL, verbose level 2 MUST be used only for debugging and
851 //it is highly suggested to use this mode only for debugging digits files
852 //reasonably small, because otherwise the size of the txt files can reach
853 //quickly several MB wasting time and disk space.
854 buffer->SetVerbose(0);
856 Int_t nEntries = Int_t(digits->GetEntries());
857 Int_t previousSector = -1;
859 for (Int_t i = 0; i < nEntries; i++) {
862 fTPCParam->AdjustSectorRow(digarr.GetID(), sector, row);
863 if(previousSector != sector) {
865 previousSector = sector;
868 if (sector < 36) { //inner sector [0;35]
870 //the whole row is written into the output file
871 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
872 sector, subSector, row);
874 //only the pads in the range [37;48] are written into the output file
875 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 1,
876 sector, subSector, row);
878 //only the pads outside the range [37;48] are written into the output file
879 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 2,
880 sector, subSector, row);
883 } else { //outer sector [36;71]
884 if (row == 54) subSector = 2;
885 if ((row != 27) && (row != 76)) {
886 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
887 sector, subSector, row);
888 } else if (row == 27) {
889 //only the pads outside the range [43;46] are written into the output file
890 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 2,
891 sector, subSector, row);
893 //only the pads in the range [43;46] are written into the output file
894 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 1,
895 sector, subSector, row);
896 } else if (row == 76) {
897 //only the pads outside the range [33;88] are written into the output file
898 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 2,
899 sector, subSector, row);
901 //only the pads in the range [33;88] are written into the output file
902 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 1,
903 sector, subSector, row);
909 fLoader->UnloadDigits();
911 AliTPCDDLRawData rawWriter;
912 rawWriter.SetVerbose(0);
914 rawWriter.RawData(fileName);
915 gSystem->Unlink(fileName);
920 //_____________________________________________________________________________
921 Bool_t AliTPC::Raw2SDigits(AliRawReader* rawReader){
922 // Converts the TPC raw data into summable digits
923 // The method is used for merging simulated and
925 if (fLoader->TreeS() == 0x0 ) {
926 fLoader->MakeTree("S");
929 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
931 //setup TPCDigitsArray
932 if(GetDigitsArray()) delete GetDigitsArray();
934 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
935 arr->SetClass("AliSimDigits");
936 arr->Setup(fTPCParam);
937 arr->MakeTree(fLoader->TreeS());
941 // set zero suppression to "0"
942 fTPCParam->SetZeroSup(0);
945 const Int_t kmaxTime = fTPCParam->GetMaxTBin();
946 const Int_t kNIS = fTPCParam->GetNInnerSector();
947 const Int_t kNOS = fTPCParam->GetNOuterSector();
948 const Int_t kNS = kNIS + kNOS;
951 AliTPCROC * roc = AliTPCROC::Instance();
952 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
953 Int_t nPadsMax = roc->GetNPads(roc->GetNSector()-1,nRowsMax-1);
954 Short_t** allBins = new Short_t*[nRowsMax];
955 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
956 Int_t maxBin = kmaxTime*nPadsMax;
957 allBins[iRow] = new Short_t[maxBin];
958 memset(allBins[iRow],0,sizeof(Short_t)*maxBin);
961 for(Int_t iSector = 0; iSector < kNS; iSector++) {
963 Int_t nRows = fTPCParam->GetNRow(iSector);
964 Int_t nDDLs = 0, indexDDL = 0;
965 if (iSector < kNIS) {
967 indexDDL = iSector * 2;
971 indexDDL = (iSector-kNIS) * 4 + kNIS * 2;
974 // Load the raw data for corresponding DDLs
977 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();
978 AliTPCRawStreamV3 input(rawReader,(AliAltroMapping**)mapping);
979 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
982 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
983 Int_t maxBin = kmaxTime*nPadsMax;
984 memset(allBins[iRow],0,sizeof(Short_t)*maxBin);
987 // Begin loop over altro data
988 while (input.NextDDL()) {
990 if (input.GetSector() != iSector)
991 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",iSector,input.GetSector()));
994 while ( input.NextChannel() ) {
996 Int_t iRow = input.GetRow();
997 if (iRow < 0 || iRow >= nRows)
998 AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
1000 Int_t iPad = input.GetPad();
1002 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
1004 if (iPad < 0 || iPad >= maxPad)
1005 AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
1006 iPad, 0, maxPad -1));
1009 while ( input.NextBunch() ){
1010 Int_t startTbin = (Int_t)input.GetStartTimeBin();
1011 Int_t bunchlength = (Int_t)input.GetBunchLength();
1012 const UShort_t *sig = input.GetSignals();
1013 for (Int_t iTime = 0; iTime<bunchlength; iTime++){
1014 Int_t iTimeBin=startTbin-iTime;
1015 if ( iTimeBin < 0 || iTimeBin >= kmaxTime) {
1017 //AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
1018 // iTimeBin, 0, kmaxTime -1));
1021 Int_t maxBin = kmaxTime*maxPad;
1022 if (((iPad*kmaxTime+iTimeBin) >= maxBin) ||
1023 ((iPad*kmaxTime+iTimeBin) < 0))
1024 AliFatal(Form("Index outside the allowed range"
1025 " Sector=%d Row=%d Pad=%d Timebin=%d"
1026 " (Max.index=%d)",iSector,iRow,iPad,iTimeBin,maxBin));
1027 allBins[iRow][iPad*kmaxTime+iTimeBin] = sig[iTime];
1030 } // End loop over altro data
1033 // Now fill the digits array
1034 if (fDigitsArray->GetTree()==0) {
1035 AliFatal("Tree not set in fDigitsArray");
1038 for (Int_t iRow = 0; iRow < nRows; iRow++) {
1039 AliDigits * dig = fDigitsArray->CreateRow(iSector,iRow);
1040 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
1041 for(Int_t iPad = 0; iPad < maxPad; iPad++) {
1042 for(Int_t iTimeBin = 0; iTimeBin < kmaxTime; iTimeBin++) {
1043 Short_t q = allBins[iRow][iPad*kmaxTime + iTimeBin];
1044 if (q <= 0) continue;
1046 dig->SetDigitFast((Short_t)q,iTimeBin,iPad);
1049 fDigitsArray->StoreRow(iSector,iRow);
1050 Int_t ndig = dig->GetDigitSize();
1053 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1054 iSector,iRow,ndig));
1056 fDigitsArray->ClearRow(iSector,iRow);
1058 } // end of the sector digitization
1060 // get LHC clock phase from the digits tree
1062 TParameter<float> *ph;
1064 TTree *digtree = fLoader->TreeD();
1066 if(digtree){ // if TreeD exists
1067 ph = (TParameter<float>*)digtree->GetUserInfo()->FindObject("lhcphase0");
1068 phase = ph->GetVal();
1070 else{ //TreeD does not exist
1074 // store lhc clock phase in S-digits tree
1076 fLoader->TreeS()->GetUserInfo()->Add(new TParameter<float>("lhcphase0",phase));
1078 fLoader->WriteSDigits("OVERWRITE");
1080 if(GetDigitsArray()) delete GetDigitsArray();
1081 SetDigitsArray(0x0);
1084 for (Int_t iRow = 0; iRow < nRowsMax; iRow++)
1085 delete [] allBins[iRow];
1091 //______________________________________________________________________
1092 AliDigitizer* AliTPC::CreateDigitizer(AliRunDigitizer* manager) const
1094 return new AliTPCDigitizer(manager);
1097 void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
1099 //create digits from summable digits
1100 GenerNoise(500000); //create teble with noise
1102 //conect tree with sSDigits
1103 TTree *t = fLoader->TreeS();
1106 fLoader->LoadSDigits("READ");
1107 t = fLoader->TreeS();
1109 AliError("Can not get input TreeS");
1114 if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
1116 AliSimDigits digarr, *dummy=&digarr;
1117 TBranch* sdb = t->GetBranch("Segment");
1119 AliError("Can not find branch with segments in TreeS.");
1123 sdb->SetAddress(&dummy);
1125 Stat_t nentries = t->GetEntries();
1127 // set zero suppression
1129 fTPCParam->SetZeroSup(2);
1131 // get zero suppression
1133 Int_t zerosup = fTPCParam->GetZeroSup();
1135 //make tree with digits
1137 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1138 arr->SetClass("AliSimDigits");
1139 arr->Setup(fTPCParam);
1140 arr->MakeTree(fLoader->TreeD());
1142 AliTPCParam * par = fTPCParam;
1144 //Loop over segments of the TPC
1146 for (Int_t n=0; n<nentries; n++) {
1149 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1150 AliWarning(Form("Invalid segment ID ! %d",digarr.GetID()));
1153 if (!IsSectorActive(sec)) continue;
1155 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1156 Int_t nrows = digrow->GetNRows();
1157 Int_t ncols = digrow->GetNCols();
1159 digrow->ExpandBuffer();
1160 digarr.ExpandBuffer();
1161 digrow->ExpandTrackBuffer();
1162 digarr.ExpandTrackBuffer();
1165 Short_t * pamp0 = digarr.GetDigits();
1166 Int_t * ptracks0 = digarr.GetTracks();
1167 Short_t * pamp1 = digrow->GetDigits();
1168 Int_t * ptracks1 = digrow->GetTracks();
1169 Int_t nelems =nrows*ncols;
1170 Int_t saturation = fTPCParam->GetADCSat() - 1;
1171 //use internal structure of the AliDigits - for speed reason
1172 //if you cahnge implementation
1173 //of the Alidigits - it must be rewriten -
1174 for (Int_t i= 0; i<nelems; i++){
1175 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1177 if (q>saturation) q=saturation;
1180 ptracks1[0]=ptracks0[0];
1181 ptracks1[nelems]=ptracks0[nelems];
1182 ptracks1[2*nelems]=ptracks0[2*nelems];
1190 arr->StoreRow(sec,row);
1191 arr->ClearRow(sec,row);
1196 fLoader->WriteDigits("OVERWRITE");
1200 //__________________________________________________________________
1201 void AliTPC::SetDefaults(){
1203 // setting the defaults
1206 // Set response functions
1209 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1211 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1213 // AliInfo("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...");
1215 // param = new AliTPCParamSR();
1218 // param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1220 param = (AliTPCParamSR*)AliTPCcalibDB::Instance()->GetParameters();
1221 if (!param->IsGeoRead()){
1223 // read transformation matrices for gGeoManager
1225 param->ReadGeoMatrices();
1228 AliFatal("No TPC parameters found");
1232 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1233 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1234 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1237 //AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1238 //rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1239 //rf->SetOffset(3*param->GetZSigma());
1244 char strgamma4[1000];
1245 sprintf(strgamma4,"AliTPCRF1D::Gamma4((x-0.135+%f)*%f,55,160)",3*param->GetZSigma(), 1000000000*param->GetTSample()/param->GetZWidth());
1247 TF1 * fgamma4 = new TF1("fgamma4",strgamma4, -1,1);
1248 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE,1000);
1249 rf->SetParam(fgamma4,param->GetZWidth(), 1,0.2);
1250 rf->SetOffset(3*param->GetZSigma());
1253 TDirectory *savedir=gDirectory;
1254 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1256 AliFatal("Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !");
1259 prfinner->Read("prf_07504_Gati_056068_d02");
1260 //PH Set different names
1261 s=prfinner->GetGRF()->GetName();
1263 prfinner->GetGRF()->SetName(s.Data());
1265 prfouter1->Read("prf_10006_Gati_047051_d03");
1266 s=prfouter1->GetGRF()->GetName();
1268 prfouter1->GetGRF()->SetName(s.Data());
1270 prfouter2->Read("prf_15006_Gati_047051_d03");
1271 s=prfouter2->GetGRF()->GetName();
1273 prfouter2->GetGRF()->SetName(s.Data());
1278 param->SetInnerPRF(prfinner);
1279 param->SetOuter1PRF(prfouter1);
1280 param->SetOuter2PRF(prfouter2);
1281 param->SetTimeRF(rf);
1291 //__________________________________________________________________
1292 void AliTPC::Hits2Digits()
1295 // creates digits from hits
1297 if (!fTPCParam->IsGeoRead()){
1299 // read transformation matrices for gGeoManager
1301 fTPCParam->ReadGeoMatrices();
1304 fLoader->LoadHits("read");
1305 fLoader->LoadDigits("recreate");
1306 AliRunLoader* runLoader = fLoader->GetRunLoader();
1308 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1309 //PH runLoader->GetEvent(iEvent);
1310 Hits2Digits(iEvent);
1313 fLoader->UnloadHits();
1314 fLoader->UnloadDigits();
1316 //__________________________________________________________________
1317 void AliTPC::Hits2Digits(Int_t eventnumber)
1319 //----------------------------------------------------
1320 // Loop over all sectors for a single event
1321 //----------------------------------------------------
1322 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1323 rl->GetEvent(eventnumber);
1325 if (fLoader->TreeH() == 0x0) {
1326 if(fLoader->LoadHits()) {
1327 AliError("Can not load hits.");
1332 if (fLoader->TreeD() == 0x0 ) {
1333 fLoader->MakeTree("D");
1334 if (fLoader->TreeD() == 0x0 ) {
1335 AliError("Can not get TreeD");
1340 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1341 GenerNoise(500000); //create teble with noise
1343 //setup TPCDigitsArray
1345 if(GetDigitsArray()) delete GetDigitsArray();
1347 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1348 arr->SetClass("AliSimDigits");
1349 arr->Setup(fTPCParam);
1351 arr->MakeTree(fLoader->TreeD());
1352 SetDigitsArray(arr);
1354 fDigitsSwitch=0; // standard digits
1355 // here LHC clock phase
1357 switch (fLHCclockPhaseSw){
1364 lhcph = (Int_t)(gRandom->Rndm()/0.25);
1368 // not implemented yet
1371 // adding phase to the TreeD user info
1372 fLoader->TreeD()->GetUserInfo()->Add(new TParameter<float>("lhcphase0",lhcph));
1374 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1375 if (IsSectorActive(isec)) {
1376 AliDebug(1,Form("Hits2Digits: Sector %d is active.",isec));
1377 Hits2DigitsSector(isec);
1380 AliDebug(1,Form("Hits2Digits: Sector %d is NOT active.",isec));
1383 fLoader->WriteDigits("OVERWRITE");
1385 //this line prevents the crash in the similar one
1386 //on the beginning of this method
1387 //destructor attempts to reset the tree, which is deleted by the loader
1388 //need to be redesign
1389 if(GetDigitsArray()) delete GetDigitsArray();
1390 SetDigitsArray(0x0);
1394 //__________________________________________________________________
1395 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1398 //-----------------------------------------------------------
1399 // summable digits - 16 bit "ADC", no noise, no saturation
1400 //-----------------------------------------------------------
1402 //----------------------------------------------------
1403 // Loop over all sectors for a single event
1404 //----------------------------------------------------
1406 AliRunLoader* rl = fLoader->GetRunLoader();
1408 rl->GetEvent(eventnumber);
1409 if (fLoader->TreeH() == 0x0) {
1410 if(fLoader->LoadHits()) {
1411 AliError("Can not load hits.");
1418 if (fLoader->TreeS() == 0x0 ) {
1419 fLoader->MakeTree("S");
1422 if(fDefaults == 0) SetDefaults();
1424 GenerNoise(500000); //create table with noise
1425 //setup TPCDigitsArray
1427 if(GetDigitsArray()) delete GetDigitsArray();
1430 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1431 arr->SetClass("AliSimDigits");
1432 arr->Setup(fTPCParam);
1433 arr->MakeTree(fLoader->TreeS());
1435 SetDigitsArray(arr);
1437 fDigitsSwitch=1; // summable digits
1439 // set zero suppression to "0"
1440 // here LHC clock phase
1442 switch (fLHCclockPhaseSw){
1449 lhcph = (Int_t)(gRandom->Rndm()/0.25);
1453 // not implemented yet
1456 // adding phase to the TreeS user info
1458 fLoader->TreeS()->GetUserInfo()->Add(new TParameter<float>("lhcphase0",lhcph));
1460 fTPCParam->SetZeroSup(0);
1462 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1463 if (IsSectorActive(isec)) {
1464 Hits2DigitsSector(isec);
1467 fLoader->WriteSDigits("OVERWRITE");
1469 //this line prevents the crash in the similar one
1470 //on the beginning of this method
1471 //destructor attempts to reset the tree, which is deleted by the loader
1472 //need to be redesign
1473 if(GetDigitsArray()) delete GetDigitsArray();
1474 SetDigitsArray(0x0);
1476 //__________________________________________________________________
1478 void AliTPC::Hits2SDigits()
1481 //-----------------------------------------------------------
1482 // summable digits - 16 bit "ADC", no noise, no saturation
1483 //-----------------------------------------------------------
1485 if (!fTPCParam->IsGeoRead()){
1487 // read transformation matrices for gGeoManager
1489 fTPCParam->ReadGeoMatrices();
1492 fLoader->LoadHits("read");
1493 fLoader->LoadSDigits("recreate");
1494 AliRunLoader* runLoader = fLoader->GetRunLoader();
1496 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1497 runLoader->GetEvent(iEvent);
1500 Hits2SDigits2(iEvent);
1503 fLoader->UnloadHits();
1504 fLoader->UnloadSDigits();
1505 if (fDebugStreamer) {
1506 delete fDebugStreamer;
1510 //_____________________________________________________________________________
1512 void AliTPC::Hits2DigitsSector(Int_t isec)
1514 //-------------------------------------------------------------------
1515 // TPC conversion from hits to digits.
1516 //-------------------------------------------------------------------
1518 //-----------------------------------------------------------------
1519 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1520 //-----------------------------------------------------------------
1522 //-------------------------------------------------------
1523 // Get the access to the track hits
1524 //-------------------------------------------------------
1526 // check if the parameters are set - important if one calls this method
1527 // directly, not from the Hits2Digits
1529 if(fDefaults == 0) SetDefaults();
1531 TTree *tH = fLoader->TreeH(); // pointer to the hits tree
1533 AliFatal("Can not find TreeH in folder");
1537 Stat_t ntracks = tH->GetEntries();
1543 Int_t nrows =fTPCParam->GetNRow(isec);
1545 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1547 MakeSector(isec,nrows,tH,ntracks,row);
1549 //--------------------------------------------------------
1550 // Digitize this sector, row by row
1551 // row[i] is the pointer to the TObjArray of TVectors,
1552 // each one containing electrons accepted on this
1553 // row, assigned into tracks
1554 //--------------------------------------------------------
1558 if (fDigitsArray->GetTree()==0) {
1559 AliFatal("Tree not set in fDigitsArray");
1562 for (i=0;i<nrows;i++){
1564 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1566 DigitizeRow(i,isec,row);
1568 fDigitsArray->StoreRow(isec,i);
1570 Int_t ndig = dig->GetDigitSize();
1573 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1576 fDigitsArray->ClearRow(isec,i);
1579 } // end of the sector digitization
1581 for(i=0;i<nrows+2;i++){
1586 delete [] row; // delete the array of pointers to TObjArray-s
1589 } // end of Hits2DigitsSector
1592 //_____________________________________________________________________________
1593 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1595 //-----------------------------------------------------------
1596 // Single row digitization, coupling from the neighbouring
1597 // rows taken into account
1598 //-----------------------------------------------------------
1600 //-----------------------------------------------------------------
1601 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1602 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1603 //-----------------------------------------------------------------
1605 Float_t zerosup = fTPCParam->GetZeroSup();
1606 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor();
1607 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
1608 AliTPCCalROC * gainROC = gainTPC->GetCalROC(isec); // pad gains per given sector
1609 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(isec); // noise per given sector
1612 fCurrentIndex[1]= isec;
1615 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1616 Int_t nofTbins = fTPCParam->GetMaxTBin();
1617 Int_t indexRange[4];
1619 // Integrated signal for this row
1620 // and a single track signal
1623 TMatrixF *m1 = new TMatrixF(0,nofPads,0,nofTbins); // integrated
1624 TMatrixF *m2 = new TMatrixF(0,nofPads,0,nofTbins); // single
1626 TMatrixF &total = *m1;
1628 // Array of pointers to the label-signal list
1630 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1631 Float_t **pList = new Float_t* [nofDigits];
1635 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1641 for (Int_t row= row1;row<=row2;row++){
1642 Int_t nTracks= rows[row]->GetEntries();
1643 for (i1=0;i1<nTracks;i1++){
1644 fCurrentIndex[2]= row;
1645 fCurrentIndex[3]=irow+1;
1647 m2->Zero(); // clear single track signal matrix
1648 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1649 GetList(trackLabel,nofPads,m2,indexRange,pList);
1651 else GetSignal(rows[row],i1,0,m1,indexRange);
1657 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1659 Float_t fzerosup = zerosup+0.5;
1660 for(Int_t it=0;it<nofTbins;it++){
1661 for(Int_t ip=0;ip<nofPads;ip++){
1663 Float_t q=total(ip,it);
1664 if(fDigitsSwitch == 0){
1665 Float_t gain = gainROC->GetValue(irow,ip); // get gain for given - pad-row pad
1666 Float_t noisePad = noiseROC->GetValue(irow,ip);
1669 q+=GetNoise()*noisePad;
1670 if(q <=fzerosup) continue; // do not fill zeros
1672 if(q >= fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat() - 1; // saturation
1677 if(q <= 0.) continue; // do not fill zeros
1678 if(q>2000.) q=2000.;
1684 // "real" signal or electronic noise (list = -1)?
1687 for(Int_t j1=0;j1<3;j1++){
1688 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1693 <A NAME="AliDigits"></A>
1694 using of AliDigits object
1697 dig->SetDigitFast((Short_t)q,it,ip);
1698 if (fDigitsArray->IsSimulated()) {
1699 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1700 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1701 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1704 } // end of loop over time buckets
1705 } // end of lop over pads
1711 // glitch filters if normal simulated digits
1713 if(!fDigitsSwitch) ((AliSimDigits*)dig)->GlitchFilter();
1715 // This row has been digitized, delete nonused stuff
1718 for(lp=0;lp<nofDigits;lp++){
1719 if(pList[lp]) delete [] pList[lp];
1727 } // end of DigitizeRow
1729 //_____________________________________________________________________________
1731 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1732 TMatrixF *m1, TMatrixF *m2,Int_t *indexRange)
1735 //---------------------------------------------------------------
1736 // Calculates 2-D signal (pad,time) for a single track,
1737 // returns a pointer to the signal matrix and the track label
1738 // No digitization is performed at this level!!!
1739 //---------------------------------------------------------------
1741 //-----------------------------------------------------------------
1742 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1743 // Modified: Marian Ivanov
1744 //-----------------------------------------------------------------
1748 tv = (TVector*)p1->At(ntr); // pointer to a track
1751 Float_t label = v(0);
1752 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1))/2;
1754 Int_t nElectrons = (tv->GetNrows()-1)/5;
1755 indexRange[0]=9999; // min pad
1756 indexRange[1]=-1; // max pad
1757 indexRange[2]=9999; //min time
1758 indexRange[3]=-1; // max time
1760 TMatrixF &signal = *m1;
1761 TMatrixF &total = *m2;
1763 // Get LHC clock phase
1765 TParameter<float> *ph;
1766 if(fDigitsSwitch){// s-digits
1767 ph = (TParameter<float>*)fLoader->TreeS()->GetUserInfo()->FindObject("lhcphase0");
1769 else{ // normal digits
1770 ph = (TParameter<float>*)fLoader->TreeD()->GetUserInfo()->FindObject("lhcphase0");
1772 // Loop over all electrons
1774 for(Int_t nel=0; nel<nElectrons; nel++){
1776 Float_t aval = v(idx+4);
1777 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1778 Float_t xyz[4]={v(idx+1),v(idx+2),v(idx+3),v(idx+5)};
1779 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,
1780 fCurrentIndex[3],ph->GetVal());
1782 Int_t *index = fTPCParam->GetResBin(0);
1783 Float_t *weight = & (fTPCParam->GetResWeight(0));
1785 if (n>0) for (Int_t i =0; i<n; i++){
1786 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1789 Int_t time=index[2];
1790 Float_t qweight = *(weight)*eltoadcfac;
1792 if (m1!=0) signal(pad,time)+=qweight;
1793 total(pad,time)+=qweight;
1794 if (indexRange[0]>pad) indexRange[0]=pad;
1795 if (indexRange[1]<pad) indexRange[1]=pad;
1796 if (indexRange[2]>time) indexRange[2]=time;
1797 if (indexRange[3]<time) indexRange[3]=time;
1804 } // end of loop over electrons
1806 return label; // returns track label when finished
1809 //_____________________________________________________________________________
1810 void AliTPC::GetList(Float_t label,Int_t np,TMatrixF *m,
1811 Int_t *indexRange, Float_t **pList)
1813 //----------------------------------------------------------------------
1814 // Updates the list of tracks contributing to digits for a given row
1815 //----------------------------------------------------------------------
1817 //-----------------------------------------------------------------
1818 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1819 //-----------------------------------------------------------------
1821 TMatrixF &signal = *m;
1823 // lop over nonzero digits
1825 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1826 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1829 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1831 if(signal(ip,it)<0.5) continue;
1833 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1835 if(!pList[globalIndex]){
1838 // Create new list (6 elements - 3 signals and 3 labels),
1841 pList[globalIndex] = new Float_t [6];
1845 *pList[globalIndex] = -1.;
1846 *(pList[globalIndex]+1) = -1.;
1847 *(pList[globalIndex]+2) = -1.;
1848 *(pList[globalIndex]+3) = -1.;
1849 *(pList[globalIndex]+4) = -1.;
1850 *(pList[globalIndex]+5) = -1.;
1852 *pList[globalIndex] = label;
1853 *(pList[globalIndex]+3) = signal(ip,it);
1857 // check the signal magnitude
1859 Float_t highest = *(pList[globalIndex]+3);
1860 Float_t middle = *(pList[globalIndex]+4);
1861 Float_t lowest = *(pList[globalIndex]+5);
1864 // compare the new signal with already existing list
1867 if(signal(ip,it)<lowest) continue; // neglect this track
1871 if (signal(ip,it)>highest){
1872 *(pList[globalIndex]+5) = middle;
1873 *(pList[globalIndex]+4) = highest;
1874 *(pList[globalIndex]+3) = signal(ip,it);
1876 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1877 *(pList[globalIndex]+1) = *pList[globalIndex];
1878 *pList[globalIndex] = label;
1880 else if (signal(ip,it)>middle){
1881 *(pList[globalIndex]+5) = middle;
1882 *(pList[globalIndex]+4) = signal(ip,it);
1884 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1885 *(pList[globalIndex]+1) = label;
1888 *(pList[globalIndex]+5) = signal(ip,it);
1889 *(pList[globalIndex]+2) = label;
1893 } // end of loop over pads
1894 } // end of loop over time bins
1897 //___________________________________________________________________
1898 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1899 Stat_t ntracks,TObjArray **row)
1902 //-----------------------------------------------------------------
1903 // Prepares the sector digitization, creates the vectors of
1904 // tracks for each row of this sector. The track vector
1905 // contains the track label and the position of electrons.
1906 //-----------------------------------------------------------------
1909 // The trasport of the electrons through TPC drift volume
1910 // Drift (drift velocity + velocity map(not yet implemented)))
1911 // Application of the random processes (diffusion, gas gain)
1912 // Systematic effects (ExB effect in drift volume + ROCs)
1915 // Loop over primary electrons:
1916 // Creation of the secondary electrons
1917 // Loop over electrons (primary+ secondaries)
1918 // Global coordinate frame:
1919 // 1. Skip electrons if attached
1920 // 2. ExB effect in drift volume
1921 // a.) Simulation calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1922 // b.) Reconstruction - calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1923 // 3. Generation of gas gain (Random - Exponential distribution)
1924 // 4. TransportElectron function (diffusion)
1926 // 5. Conversion to the local coordinate frame pad-row, pad, timebin
1927 // 6. Apply Time0 shift - AliTPCCalPad class
1928 // a.) Plus sign in simulation
1929 // b.) Minus sign in reconstruction
1932 //-----------------------------------------------------------------
1933 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1934 // Origin: Marian Ivanov, marian.ivanov@cern.ch
1935 //-----------------------------------------------------------------
1936 AliTPCcalibDB* const calib=AliTPCcalibDB::Instance();
1937 if (gAlice){ // Set correctly the magnetic field in the ExB calculation
1938 if (!calib->GetExB()){
1939 AliMagF * field = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField());
1941 calib->SetExBField(field);
1946 Float_t gasgain = fTPCParam->GetGasGain();
1947 gasgain = gasgain/fGainFactor;
1951 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1954 if (fHitType>1) branch = TH->GetBranch("TPC2");
1955 else branch = TH->GetBranch("TPC");
1958 //----------------------------------------------
1959 // Create TObjArray-s, one for each row,
1960 // each TObjArray will store the TVectors
1961 // of electrons, one TVectors per each track.
1962 //----------------------------------------------
1964 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
1965 TVector **tracks = new TVector* [nrows+2]; //pointers to the track vectors
1967 for(i=0; i<nrows+2; i++){
1968 row[i] = new TObjArray;
1975 //--------------------------------------------------------------------
1976 // Loop over tracks, the "track" contains the full history
1977 //--------------------------------------------------------------------
1979 Int_t previousTrack,currentTrack;
1980 previousTrack = -1; // nothing to store so far!
1982 for(Int_t track=0;track<ntracks;track++){
1983 Bool_t isInSector=kTRUE;
1985 isInSector = TrackInVolume(isec,track);
1986 if (!isInSector) continue;
1988 branch->GetEntry(track); // get next track
1992 tpcHit = (AliTPChit*)FirstHit(-1);
1994 //--------------------------------------------------------------
1996 //--------------------------------------------------------------
2001 Int_t sector=tpcHit->fSector; // sector number
2003 tpcHit = (AliTPChit*) NextHit();
2007 // Remove hits which arrive before the TPC opening gate signal
2008 if(((fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
2009 /fTPCParam->GetDriftV()+tpcHit->Time())<fTPCParam->GetGateDelay()) {
2010 tpcHit = (AliTPChit*) NextHit();
2014 currentTrack = tpcHit->Track(); // track number
2016 if(currentTrack != previousTrack){
2018 // store already filled fTrack
2020 for(i=0;i<nrows+2;i++){
2021 if(previousTrack != -1){
2022 if(nofElectrons[i]>0){
2023 TVector &v = *tracks[i];
2024 v(0) = previousTrack;
2025 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
2026 row[i]->Add(tracks[i]);
2029 delete tracks[i]; // delete empty TVector
2035 tracks[i] = new TVector(601); // TVectors for the next fTrack
2037 } // end of loop over rows
2039 previousTrack=currentTrack; // update track label
2042 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2044 //---------------------------------------------------
2045 // Calculate the electron attachment probability
2046 //---------------------------------------------------
2049 Float_t time = 1.e6*(fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
2050 /fTPCParam->GetDriftV();
2052 Float_t attProb = fTPCParam->GetAttCoef()*
2053 fTPCParam->GetOxyCont()*time; // fraction!
2055 //-----------------------------------------------
2056 // Loop over electrons
2057 //-----------------------------------------------
2060 for(Int_t nel=0;nel<qI;nel++){
2061 // skip if electron lost due to the attachment
2062 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2067 Double_t dxyz0[3],dxyz1[3];
2068 dxyz0[0]=tpcHit->X();
2069 dxyz0[1]=tpcHit->Y();
2070 dxyz0[2]=tpcHit->Z();
2071 if (calib->GetExB()){
2072 calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
2074 AliError("Not valid ExB calibration");
2075 dxyz1[0]=tpcHit->X();
2076 dxyz1[1]=tpcHit->Y();
2077 dxyz1[2]=tpcHit->Z();
2085 // protection for the nonphysical avalanche size (10**6 maximum)
2087 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2088 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2091 TransportElectron(xyz,index);
2093 Int_t padrow = fTPCParam->GetPadRow(xyz,index);
2095 // Add Time0 correction due unisochronity
2096 // xyz[0] - pad row coordinate
2097 // xyz[1] - pad coordinate
2098 // xyz[2] - is in now time bin coordinate system
2099 Float_t correction =0;
2100 if (calib->GetPadTime0()){
2101 if (!calib->GetPadTime0()->GetCalROC(isec)) continue;
2102 Int_t npads = fTPCParam->GetNPads(isec,padrow);
2103 // Int_t pad = TMath::Nint(xyz[1]+fTPCParam->GetNPads(isec,TMath::Nint(xyz[0]))*0.5);
2104 // pad numbering from -npads/2 .. npads/2-1
2105 Int_t pad = TMath::Nint(xyz[1]+npads/2);
2107 if (pad>=npads) pad=npads-1;
2108 correction = calib->GetPadTime0()->GetCalROC(isec)->GetValue(padrow,pad);
2109 // printf("%d\t%d\t%d\t%f\n",isec,padrow,pad,correction);
2110 if (fDebugStreamer){
2111 (*fDebugStreamer)<<"Time0"<<
2119 "cor="<<correction<<
2124 xyz[2]+=fTPCParam->GetNTBinsL1(); // adding Level 1 time bin offset
2126 // Electron track time (for pileup simulation)
2127 xyz[2]+=tpcHit->Time()/fTPCParam->GetTSample(); // adding time of flight
2131 // row 0 - cross talk from the innermost row
2132 // row fNRow+1 cross talk from the outermost row
2133 rowNumber = index[2]+1;
2134 //transform position to local digit coordinates
2135 //relative to nearest pad row
2136 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2138 if (isec <fTPCParam->GetNInnerSector()) {
2139 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2140 y1 = fTPCParam->GetYInner(rowNumber);
2143 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2144 y1 = fTPCParam->GetYOuter(rowNumber);
2146 // gain inefficiency at the wires edges - linear
2149 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.)); */
2151 nofElectrons[rowNumber]++;
2152 //----------------------------------
2153 // Expand vector if necessary
2154 //----------------------------------
2155 if(nofElectrons[rowNumber]>120){
2156 Int_t range = tracks[rowNumber]->GetNrows();
2157 if((nofElectrons[rowNumber])>(range-1)/5){
2159 tracks[rowNumber]->ResizeTo(range+500); // Add 100 electrons
2163 TVector &v = *tracks[rowNumber];
2164 Int_t idx = 5*nofElectrons[rowNumber]-4;
2165 Real_t * position = &(((TVector&)v)(idx)); //make code faster
2166 memcpy(position,xyz,5*sizeof(Float_t));
2168 } // end of loop over electrons
2170 tpcHit = (AliTPChit*)NextHit();
2172 } // end of loop over hits
2173 } // end of loop over tracks
2176 // store remaining track (the last one) if not empty
2179 for(i=0;i<nrows+2;i++){
2180 if(nofElectrons[i]>0){
2181 TVector &v = *tracks[i];
2182 v(0) = previousTrack;
2183 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
2184 row[i]->Add(tracks[i]);
2193 delete [] nofElectrons;
2195 } // end of MakeSector
2198 //_____________________________________________________________________________
2202 // Initialise TPC detector after definition of geometry
2204 AliDebug(1,"*********************************************");
2207 //_____________________________________________________________________________
2208 void AliTPC::ResetDigits()
2211 // Reset number of digits and the digits array for this detector
2214 if (fDigits) fDigits->Clear();
2219 //_____________________________________________________________________________
2220 void AliTPC::SetSens(Int_t sens)
2223 //-------------------------------------------------------------
2224 // Activates/deactivates the sensitive strips at the center of
2225 // the pad row -- this is for the space-point resolution calculations
2226 //-------------------------------------------------------------
2228 //-----------------------------------------------------------------
2229 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2230 //-----------------------------------------------------------------
2236 void AliTPC::SetSide(Float_t side=0.)
2238 // choice of the TPC side
2243 //_____________________________________________________________________________
2245 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2248 // electron transport taking into account:
2250 // 2.ExB at the wires
2251 // 3. nonisochronity
2253 // xyz and index must be already transformed to system 1
2256 fTPCParam->Transform1to2(xyz,index); // mis-alignment applied in this step
2259 Float_t driftl=xyz[2];
2260 if(driftl<0.01) driftl=0.01;
2261 driftl=TMath::Sqrt(driftl);
2262 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2263 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2264 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2265 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2266 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2270 if (fTPCParam->GetMWPCReadout()==kTRUE){
2271 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2272 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2274 //add nonisochronity (not implemented yet)
2280 //______________________________________________________________________
2281 AliTPChit::AliTPChit()
2293 //_____________________________________________________________________________
2294 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits)
2295 :AliHit(shunt,track),
2302 // Creates a TPC hit object
2313 //________________________________________________________________________
2314 // Additional code because of the AliTPCTrackHitsV2
2316 void AliTPC::MakeBranch(Option_t *option)
2319 // Create a new branch in the current Root Tree
2320 // The branch of fHits is automatically split
2321 // MI change 14.09.2000
2323 if (fHitType<2) return;
2324 char branchname[10];
2325 sprintf(branchname,"%s2",GetName());
2327 // Get the pointer to the header
2328 const char *cH = strstr(option,"H");
2330 if (fTrackHits && fLoader->TreeH() && cH && fHitType&4) {
2331 AliDebug(1,"Making branch for Type 4 Hits");
2332 fLoader->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2335 // if (fTrackHitsOld && fLoader->TreeH() && cH && fHitType&2) {
2336 // AliDebug(1,"Making branch for Type 2 Hits");
2337 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2338 // fLoader->TreeH(),fBufferSize,99);
2339 // fLoader->TreeH()->GetListOfBranches()->Add(branch);
2343 void AliTPC::SetTreeAddress()
2345 //Sets tree address for hits
2347 if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2348 AliDetector::SetTreeAddress();
2350 if (fHitType>1) SetTreeAddress2();
2353 void AliTPC::SetTreeAddress2()
2356 // Set branch address for the TrackHits Tree
2361 char branchname[20];
2362 sprintf(branchname,"%s2",GetName());
2364 // Branch address for hit tree
2365 TTree *treeH = fLoader->TreeH();
2366 if ((treeH)&&(fHitType&4)) {
2367 branch = treeH->GetBranch(branchname);
2369 branch->SetAddress(&fTrackHits);
2370 AliDebug(1,"fHitType&4 Setting");
2373 AliDebug(1,"fHitType&4 Failed (can not find branch)");
2376 // if ((treeH)&&(fHitType&2)) {
2377 // branch = treeH->GetBranch(branchname);
2379 // branch->SetAddress(&fTrackHitsOld);
2380 // AliDebug(1,"fHitType&2 Setting");
2383 // AliDebug(1,"fHitType&2 Failed (can not find branch)");
2387 void AliTPC::FinishPrimary()
2389 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2390 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2394 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2397 // add hit to the list
2401 int primary = gAlice->GetMCApp()->GetPrimary(track);
2402 gAlice->GetMCApp()->Particle(primary)->SetBit(kKeepBit);
2406 gAlice->GetMCApp()->FlagTrack(track);
2408 if (fTrackHits && fHitType&4)
2409 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2410 hits[1],hits[2],(Int_t)hits[3],hits[4]);
2411 // if (fTrackHitsOld &&fHitType&2 )
2412 // fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2413 // hits[1],hits[2],(Int_t)hits[3]);
2417 void AliTPC::ResetHits()
2419 if (fHitType&1) AliDetector::ResetHits();
2420 if (fHitType>1) ResetHits2();
2423 void AliTPC::ResetHits2()
2427 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2428 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2432 AliHit* AliTPC::FirstHit(Int_t track)
2434 if (fHitType>1) return FirstHit2(track);
2435 return AliDetector::FirstHit(track);
2437 AliHit* AliTPC::NextHit()
2442 if (fHitType>1) return NextHit2();
2444 return AliDetector::NextHit();
2447 AliHit* AliTPC::FirstHit2(Int_t track)
2450 // Initialise the hit iterator
2451 // Return the address of the first hit for track
2452 // If track>=0 the track is read from disk
2453 // while if track<0 the first hit of the current
2454 // track is returned
2457 gAlice->GetMCApp()->ResetHits();
2458 fLoader->TreeH()->GetEvent(track);
2461 if (fTrackHits && fHitType&4) {
2462 fTrackHits->First();
2463 return fTrackHits->GetHit();
2465 // if (fTrackHitsOld && fHitType&2) {
2466 // fTrackHitsOld->First();
2467 // return fTrackHitsOld->GetHit();
2473 AliHit* AliTPC::NextHit2()
2476 //Return the next hit for the current track
2479 // if (fTrackHitsOld && fHitType&2) {
2480 // fTrackHitsOld->Next();
2481 // return fTrackHitsOld->GetHit();
2485 return fTrackHits->GetHit();
2491 void AliTPC::RemapTrackHitIDs(Int_t *map)
2496 if (!fTrackHits) return;
2498 // if (fTrackHitsOld && fHitType&2){
2499 // AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2500 // for (UInt_t i=0;i<arr->GetSize();i++){
2501 // AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2502 // info->fTrackID = map[info->fTrackID];
2505 // if (fTrackHitsOld && fHitType&4){
2506 if (fTrackHits && fHitType&4){
2507 TClonesArray * arr = fTrackHits->GetArray();;
2508 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2509 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2510 info->SetTrackID(map[info->GetTrackID()]);
2515 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2517 //return bool information - is track in given volume
2518 //load only part of the track information
2519 //return true if current track is in volume
2522 // if (fTrackHitsOld && fHitType&2) {
2523 // TBranch * br = fLoader->TreeH()->GetBranch("fTrackHitsInfo");
2524 // br->GetEvent(track);
2525 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2526 // for (UInt_t j=0;j<ar->GetSize();j++){
2527 // if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2531 if (fTrackHits && fHitType&4) {
2532 TBranch * br1 = fLoader->TreeH()->GetBranch("fVolumes");
2533 TBranch * br2 = fLoader->TreeH()->GetBranch("fNVolumes");
2534 br2->GetEvent(track);
2535 br1->GetEvent(track);
2536 Int_t *volumes = fTrackHits->GetVolumes();
2537 Int_t nvolumes = fTrackHits->GetNVolumes();
2538 if (!volumes && nvolumes>0) {
2539 AliWarning(Form("Problematic track\t%d\t%d",track,nvolumes));
2542 for (Int_t j=0;j<nvolumes; j++)
2543 if (volumes[j]==id) return kTRUE;
2547 TBranch * br = fLoader->TreeH()->GetBranch("fSector");
2548 br->GetEvent(track);
2549 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2550 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2558 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
2561 fLoader = new AliTPCLoader(GetName(),topfoldername);
2565 ////////////////////////////////////////////////////////////////////////
2566 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
2568 // load TPC paarmeters from a given file or create new if the object
2569 // is not found there
2570 // 12/05/2003 This method should be moved to the AliTPCLoader
2571 // and one has to decide where to store the TPC parameters
2574 sprintf(paramName,"75x40_100x60_150x60");
2575 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
2577 AliDebugClass(1,Form("TPC parameters %s found.",paramName));
2579 AliWarningClass("TPC parameters not found. Create new (they may be incorrect)");
2580 //paramTPC = new AliTPCParamSR;
2581 paramTPC = AliTPCcalibDB::Instance()->GetParameters();
2582 if (!paramTPC->IsGeoRead()){
2584 // read transformation matrices for gGeoManager
2586 paramTPC->ReadGeoMatrices();
2592 // the older version of parameters can be accessed with this code.
2593 // In some cases, we have old parameters saved in the file but
2594 // digits were created with new parameters, it can be distinguish
2595 // by the name of TPC TreeD. The code here is just for the case
2596 // we would need to compare with old data, uncomment it if needed.
2598 // char paramName[50];
2599 // sprintf(paramName,"75x40_100x60");
2600 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
2602 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2604 // sprintf(paramName,"75x40_100x60_150x60");
2605 // paramTPC=(AliTPCParam*)in->Get(paramName);
2607 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2609 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
2611 // paramTPC = new AliTPCParamSR;