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>
56 #include "AliDigits.h"
59 #include "AliRunLoader.h"
60 #include "AliSimDigits.h"
63 #include "AliTPCDigitsArray.h"
64 #include "AliTPCLoader.h"
65 #include "AliTPCPRF2D.h"
66 #include "AliTPCParamSR.h"
67 #include "AliTPCRF1D.h"
68 #include "AliTPCTrackHitsV2.h"
69 #include "AliTrackReference.h"
72 #include "AliTPCDigitizer.h"
73 #include "AliTPCBuffer.h"
74 #include "AliTPCDDLRawData.h"
76 #include "AliTPCcalibDB.h"
77 #include "AliTPCCalPad.h"
78 #include "AliTPCCalROC.h"
79 #include "AliTPCExB.h"
80 #include "AliRawReader.h"
81 #include "AliTPCRawStream.h"
82 #include "TTreeStream.h"
85 //_____________________________________________________________________________
86 AliTPC::AliTPC():AliDetector(),
96 fPrimaryIonisation(0),
106 // Default constructor
110 // fTrackHitsOld = 0;
111 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
112 fHitType = 4; // ROOT containers
114 fHitType = 2; //default CONTAINERS - based on ROOT structure
118 //_____________________________________________________________________________
119 AliTPC::AliTPC(const char *name, const char *title)
120 : AliDetector(name,title),
130 fPrimaryIonisation(0),
140 // Standard constructor
144 // Initialise arrays of hits and digits
145 fHits = new TClonesArray("AliTPChit", 176);
146 gAlice->GetMCApp()->AddHitList(fHits);
148 fTrackHits = new AliTPCTrackHitsV2;
149 fTrackHits->SetHitPrecision(0.002);
150 fTrackHits->SetStepPrecision(0.003);
151 fTrackHits->SetMaxDistance(100);
153 //fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
154 //fTrackHitsOld->SetHitPrecision(0.002);
155 //fTrackHitsOld->SetStepPrecision(0.003);
156 //fTrackHitsOld->SetMaxDistance(100);
159 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
160 fHitType = 4; // ROOT containers
170 // Initialise color attributes
171 //PH SetMarkerColor(kYellow);
174 // Set TPC parameters
178 if (!strcmp(title,"Default")) {
179 //fTPCParam = new AliTPCParamSR;
180 fTPCParam = AliTPCcalibDB::Instance()->GetParameters();
182 AliWarning("In Config.C you must set non-default parameters.");
188 //_____________________________________________________________________________
199 delete fTrackHits; //MI 15.09.2000
200 // delete fTrackHitsOld; //MI 10.12.2001
203 delete [] fNoiseTable;
204 delete [] fActiveSectors;
205 if (fDebugStreamer) delete fDebugStreamer;
208 //_____________________________________________________________________________
209 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
212 // Add a hit to the list
215 TClonesArray &lhits = *fHits;
216 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
219 AddHit2(track,vol,hits);
222 //_____________________________________________________________________________
223 void AliTPC::CreateMaterials()
225 //-----------------------------------------------
226 // Create Materials for for TPC simulations
227 //-----------------------------------------------
229 //-----------------------------------------------------------------
230 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
231 //-----------------------------------------------------------------
233 Int_t iSXFLD=((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ();
234 Float_t sXMGMX=((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();
236 Float_t amat[5]; // atomic numbers
237 Float_t zmat[5]; // z
238 Float_t wmat[5]; // proportions
244 //***************** Gases *************************
247 //--------------------------------------------------------------
248 // gases - air and CO2
249 //--------------------------------------------------------------
265 AliMixture(10,"CO2",amat,zmat,density,2,wmat);
280 AliMixture(11,"Air",amat,zmat,density,2,wmat);
282 //----------------------------------------------------------------
284 //----------------------------------------------------------------
287 // Drift gases 1 - nonsensitive, 2 - sensitive
288 // Ne-CO2-N (85-10-5)
307 AliMixture(12,"Ne-CO2-N-1",amat,zmat,density,4,wmat);
308 AliMixture(13,"Ne-CO2-N-2",amat,zmat,density,4,wmat);
309 AliMixture(30,"Ne-CO2-N-3",amat,zmat,density,4,wmat);
310 //----------------------------------------------------------------------
312 //----------------------------------------------------------------------
334 AliMixture(14,"Kevlar",amat,zmat,density,-4,wmat);
355 AliMixture(15,"NOMEX",amat,zmat,density,-4,wmat);
373 AliMixture(16,"Makrolon",amat,zmat,density,-3,wmat);
391 AliMixture(17, "Tedlar",amat,zmat,density,-3,wmat);
409 AliMixture(18, "Mylar",amat,zmat,density,-3,wmat);
410 // material for "prepregs"
411 // Epoxy - C14 H20 O3
414 // prepreg1 60% C-fiber, 40% epoxy (vol)
429 AliMixture(19, "Prepreg1",amat,zmat,density,3,wmat);
431 //prepreg2 60% glass-fiber, 40% epoxy
450 AliMixture(20, "Prepreg2",amat,zmat,density,4,wmat);
452 //prepreg3 50% glass-fiber, 50% epoxy
471 AliMixture(21, "Prepreg3",amat,zmat,density,4,wmat);
473 // G10 60% SiO2 40% epoxy
492 AliMixture(22, "G10",amat,zmat,density,4,wmat);
501 AliMaterial(23,"Al",amat[0],zmat[0],density,999.,999.);
503 // Si (for electronics
510 AliMaterial(24,"Si",amat[0],zmat[0],density,999.,999.);
519 AliMaterial(25,"Cu",amat[0],zmat[0],density,999.,999.);
521 // Epoxy - C14 H20 O3
537 AliMixture(26,"Epoxy",amat,zmat,density,-3,wmat);
555 AliMixture(27,"Plexiglas",amat,zmat,density,-3,wmat);
563 AliMaterial(28,"C",amat[0],zmat[0],density,999.,999.);
565 // Fe (steel for the inner heat screen)
573 AliMaterial(29,"Fe",amat[0],zmat[0],density,999.,999.);
575 // Peek - (C6H4-O-OC6H4-O-C6H4-CO)n
590 AliMixture(30,"Peek",amat,zmat,density,-3,wmat);
605 AliMixture(31,"Alumina",amat,zmat,density,-2,wmat);
624 AliMixture(32,"Water",amat,zmat,density,-2,wmat);
626 //----------------------------------------------------------
627 // tracking media for gases
628 //----------------------------------------------------------
630 AliMedium(0, "Air", 11, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
631 AliMedium(1, "Ne-CO2-N-1", 12, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
632 AliMedium(2, "Ne-CO2-N-2", 13, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
633 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
634 AliMedium(20, "Ne-CO2-N-3", 30, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
635 //-----------------------------------------------------------
636 // tracking media for solids
637 //-----------------------------------------------------------
639 AliMedium(4,"Al",23,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
640 AliMedium(5,"Kevlar",14,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
641 AliMedium(6,"Nomex",15,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
642 AliMedium(7,"Makrolon",16,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
643 AliMedium(8,"Mylar",18,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
644 AliMedium(9,"Tedlar",17,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
646 AliMedium(10,"Prepreg1",19,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
647 AliMedium(11,"Prepreg2",20,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
648 AliMedium(12,"Prepreg3",21,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
649 AliMedium(13,"Epoxy",26,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
651 AliMedium(14,"Cu",25,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
652 AliMedium(15,"Si",24,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
653 AliMedium(16,"G10",22,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
654 AliMedium(17,"Plexiglas",27,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
655 AliMedium(18,"Steel",29,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
656 AliMedium(19,"Peek",30,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
657 AliMedium(21,"Alumina",31,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
658 AliMedium(22,"Water",32,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
661 void AliTPC::GenerNoise(Int_t tablesize)
664 //Generate table with noise
671 if (fNoiseTable) delete[] fNoiseTable;
672 fNoiseTable = new Float_t[tablesize];
673 fNoiseDepth = tablesize;
674 fCurrentNoise =0; //!index of the noise in the noise table
676 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
677 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
680 Float_t AliTPC::GetNoise()
682 // get noise from table
683 // if ((fCurrentNoise%10)==0)
684 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
685 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
686 return fNoiseTable[fCurrentNoise++];
687 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
691 Bool_t AliTPC::IsSectorActive(Int_t sec) const
694 // check if the sector is active
695 if (!fActiveSectors) return kTRUE;
696 else return fActiveSectors[sec];
699 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
701 // activate interesting sectors
702 SetTreeAddress();//just for security
703 if (!fActiveSectors) fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
704 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
705 for (Int_t i=0;i<n;i++)
706 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
710 void AliTPC::SetActiveSectors(Int_t flag)
713 // activate sectors which were hitted by tracks
715 SetTreeAddress();//just for security
716 if (fHitType==0) return; // if Clones hit - not short volume ID information
717 if (!fActiveSectors) fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
719 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
722 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
724 if (fLoader->TreeH() == 0x0)
726 AliFatal("Can not find TreeH in folder");
729 if (fHitType>1) branch = fLoader->TreeH()->GetBranch("TPC2");
730 else branch = fLoader->TreeH()->GetBranch("TPC");
731 Stat_t ntracks = fLoader->TreeH()->GetEntries();
732 // loop over all hits
733 AliDebug(1,Form("Got %d tracks",ntracks));
735 for(Int_t track=0;track<ntracks;track++) {
738 if (fTrackHits && fHitType&4) {
739 TBranch * br1 = fLoader->TreeH()->GetBranch("fVolumes");
740 TBranch * br2 = fLoader->TreeH()->GetBranch("fNVolumes");
741 br1->GetEvent(track);
742 br2->GetEvent(track);
743 Int_t *volumes = fTrackHits->GetVolumes();
744 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++) {
745 if (volumes[j]>-1 && volumes[j]<fTPCParam->GetNSector()) {
746 fActiveSectors[volumes[j]]=kTRUE;
749 AliError(Form("Volume %d -> sector number %d is outside (0..%d)",
752 fTPCParam->GetNSector()));
758 // if (fTrackHitsOld && fHitType&2) {
759 // TBranch * br = fLoader->TreeH()->GetBranch("fTrackHitsInfo");
760 // br->GetEvent(track);
761 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
762 // for (UInt_t j=0;j<ar->GetSize();j++){
763 // fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
772 //_____________________________________________________________________________
773 void AliTPC::Digits2Raw()
775 // convert digits of the current event to raw data
777 static const Int_t kThreshold = 0;
779 fLoader->LoadDigits();
780 TTree* digits = fLoader->TreeD();
782 AliError("No digits tree");
787 AliSimDigits* digrow = &digarr;
788 digits->GetBranch("Segment")->SetAddress(&digrow);
790 const char* fileName = "AliTPCDDL.dat";
791 AliTPCBuffer* buffer = new AliTPCBuffer(fileName);
795 // 2: txt files with digits
796 //BE CAREFUL, verbose level 2 MUST be used only for debugging and
797 //it is highly suggested to use this mode only for debugging digits files
798 //reasonably small, because otherwise the size of the txt files can reach
799 //quickly several MB wasting time and disk space.
800 buffer->SetVerbose(0);
802 Int_t nEntries = Int_t(digits->GetEntries());
803 Int_t previousSector = -1;
805 for (Int_t i = 0; i < nEntries; i++) {
808 fTPCParam->AdjustSectorRow(digarr.GetID(), sector, row);
809 if(previousSector != sector) {
811 previousSector = sector;
814 if (sector < 36) { //inner sector [0;35]
816 //the whole row is written into the output file
817 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
818 sector, subSector, row);
820 //only the pads in the range [37;48] are written into the output file
821 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 1,
822 sector, subSector, row);
824 //only the pads outside the range [37;48] are written into the output file
825 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 2,
826 sector, subSector, row);
829 } else { //outer sector [36;71]
830 if (row == 54) subSector = 2;
831 if ((row != 27) && (row != 76)) {
832 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
833 sector, subSector, row);
834 } else if (row == 27) {
835 //only the pads outside the range [43;46] are written into the output file
836 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 2,
837 sector, subSector, row);
839 //only the pads in the range [43;46] are written into the output file
840 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 1,
841 sector, subSector, row);
842 } else if (row == 76) {
843 //only the pads outside the range [33;88] are written into the output file
844 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 2,
845 sector, subSector, row);
847 //only the pads in the range [33;88] are written into the output file
848 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 1,
849 sector, subSector, row);
855 fLoader->UnloadDigits();
857 AliTPCDDLRawData rawWriter;
858 rawWriter.SetVerbose(0);
860 rawWriter.RawData(fileName);
861 gSystem->Unlink(fileName);
866 //_____________________________________________________________________________
867 Bool_t AliTPC::Raw2SDigits(AliRawReader* rawReader){
868 // Converts the TPC raw data into summable digits
869 // The method is used for merging simulated and
871 if (fLoader->TreeS() == 0x0 ) {
872 fLoader->MakeTree("S");
875 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
877 //setup TPCDigitsArray
878 if(GetDigitsArray()) delete GetDigitsArray();
880 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
881 arr->SetClass("AliSimDigits");
882 arr->Setup(fTPCParam);
883 arr->MakeTree(fLoader->TreeS());
887 // set zero suppression to "0"
888 fTPCParam->SetZeroSup(0);
891 const Int_t kmaxTime = fTPCParam->GetMaxTBin();
892 const Int_t kNIS = fTPCParam->GetNInnerSector();
893 const Int_t kNOS = fTPCParam->GetNOuterSector();
894 const Int_t kNS = kNIS + kNOS;
896 Short_t** allBins = NULL; //array which contains the data for one sector
898 for(Int_t iSector = 0; iSector < kNS; iSector++) {
900 Int_t nRows = fTPCParam->GetNRow(iSector);
901 Int_t nDDLs = 0, indexDDL = 0;
902 if (iSector < kNIS) {
904 indexDDL = iSector * 2;
908 indexDDL = (iSector-kNIS) * 4 + kNIS * 2;
911 // Loas the raw data for corresponding DDLs
913 AliTPCRawStream input(rawReader);
914 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
916 // Alocate and init the array with the sector data
917 allBins = new Short_t*[nRows];
918 for (Int_t iRow = 0; iRow < nRows; iRow++) {
919 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
920 Int_t maxBin = kmaxTime*maxPad;
921 allBins[iRow] = new Short_t[maxBin];
922 memset(allBins[iRow],0,sizeof(Short_t)*maxBin);
925 // Begin loop over altro data
926 while (input.Next()) {
928 if (input.GetSector() != iSector)
929 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",iSector,input.GetSector()));
931 Int_t iRow = input.GetRow();
932 if (iRow < 0 || iRow >= nRows)
933 AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
935 Int_t iPad = input.GetPad();
937 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
939 if (iPad < 0 || iPad >= maxPad)
940 AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
941 iPad, 0, maxPad -1));
943 Int_t iTimeBin = input.GetTime();
944 if ( iTimeBin < 0 || iTimeBin >= kmaxTime)
945 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
946 iTimeBin, 0, kmaxTime -1));
948 Int_t maxBin = kmaxTime*maxPad;
950 if (((iPad*kmaxTime+iTimeBin) >= maxBin) ||
951 ((iPad*kmaxTime+iTimeBin) < 0))
952 AliFatal(Form("Index outside the allowed range"
953 " Sector=%d Row=%d Pad=%d Timebin=%d"
954 " (Max.index=%d)",iSector,iRow,iPad,iTimeBin,maxBin));
956 allBins[iRow][iPad*kmaxTime+iTimeBin] = input.GetSignal();
958 } // End loop over altro data
960 // Now fill the digits array
961 if (fDigitsArray->GetTree()==0) {
962 AliFatal("Tree not set in fDigitsArray");
965 for (Int_t iRow = 0; iRow < nRows; iRow++) {
966 AliDigits * dig = fDigitsArray->CreateRow(iSector,iRow);
968 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
969 for(Int_t iPad = 0; iPad < maxPad; iPad++) {
970 for(Int_t iTimeBin = 0; iTimeBin < kmaxTime; iTimeBin++) {
971 Short_t q = allBins[iRow][iPad*kmaxTime + iTimeBin];
972 if (q <= 0) continue;
974 dig->SetDigitFast((Short_t)q,iTimeBin,iPad);
977 fDigitsArray->StoreRow(iSector,iRow);
978 Int_t ndig = dig->GetDigitSize();
981 Form("*** Sector, row, compressed digits %d %d %d ***\n",
984 fDigitsArray->ClearRow(iSector,iRow);
986 } // end of the sector digitization
988 for (Int_t iRow = 0; iRow < nRows; iRow++)
989 delete [] allBins[iRow];
995 fLoader->WriteSDigits("OVERWRITE");
997 if(GetDigitsArray()) delete GetDigitsArray();
1003 //______________________________________________________________________
1004 AliDigitizer* AliTPC::CreateDigitizer(AliRunDigitizer* manager) const
1006 return new AliTPCDigitizer(manager);
1009 void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
1011 //create digits from summable digits
1012 GenerNoise(500000); //create teble with noise
1014 //conect tree with sSDigits
1015 TTree *t = fLoader->TreeS();
1018 fLoader->LoadSDigits("READ");
1019 t = fLoader->TreeS();
1021 AliError("Can not get input TreeS");
1026 if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
1028 AliSimDigits digarr, *dummy=&digarr;
1029 TBranch* sdb = t->GetBranch("Segment");
1031 AliError("Can not find branch with segments in TreeS.");
1035 sdb->SetAddress(&dummy);
1037 Stat_t nentries = t->GetEntries();
1039 // set zero suppression
1041 fTPCParam->SetZeroSup(2);
1043 // get zero suppression
1045 Int_t zerosup = fTPCParam->GetZeroSup();
1047 //make tree with digits
1049 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1050 arr->SetClass("AliSimDigits");
1051 arr->Setup(fTPCParam);
1052 arr->MakeTree(fLoader->TreeD());
1054 AliTPCParam * par = fTPCParam;
1056 //Loop over segments of the TPC
1058 for (Int_t n=0; n<nentries; n++) {
1061 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1062 AliWarning(Form("Invalid segment ID ! %d",digarr.GetID()));
1065 if (!IsSectorActive(sec)) continue;
1067 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1068 Int_t nrows = digrow->GetNRows();
1069 Int_t ncols = digrow->GetNCols();
1071 digrow->ExpandBuffer();
1072 digarr.ExpandBuffer();
1073 digrow->ExpandTrackBuffer();
1074 digarr.ExpandTrackBuffer();
1077 Short_t * pamp0 = digarr.GetDigits();
1078 Int_t * ptracks0 = digarr.GetTracks();
1079 Short_t * pamp1 = digrow->GetDigits();
1080 Int_t * ptracks1 = digrow->GetTracks();
1081 Int_t nelems =nrows*ncols;
1082 Int_t saturation = fTPCParam->GetADCSat() - 1;
1083 //use internal structure of the AliDigits - for speed reason
1084 //if you cahnge implementation
1085 //of the Alidigits - it must be rewriten -
1086 for (Int_t i= 0; i<nelems; i++){
1087 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1089 if (q>saturation) q=saturation;
1092 ptracks1[0]=ptracks0[0];
1093 ptracks1[nelems]=ptracks0[nelems];
1094 ptracks1[2*nelems]=ptracks0[2*nelems];
1102 arr->StoreRow(sec,row);
1103 arr->ClearRow(sec,row);
1108 fLoader->WriteDigits("OVERWRITE");
1112 //__________________________________________________________________
1113 void AliTPC::SetDefaults(){
1115 // setting the defaults
1118 // Set response functions
1121 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1123 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1125 // AliInfo("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...");
1127 // param = new AliTPCParamSR();
1130 // param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1132 param = (AliTPCParamSR*)AliTPCcalibDB::Instance()->GetParameters();
1133 if (!param->IsGeoRead()){
1135 // read transformation matrices for gGeoManager
1137 param->ReadGeoMatrices();
1140 AliFatal("No TPC parameters found");
1144 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1145 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1146 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1149 //AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1150 //rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1151 //rf->SetOffset(3*param->GetZSigma());
1156 char strgamma4[1000];
1157 sprintf(strgamma4,"AliTPCRF1D::Gamma4((x-0.135+%f)*%f,55,160)",3*param->GetZSigma(), 1000000000*param->GetTSample()/param->GetZWidth());
1159 TF1 * fgamma4 = new TF1("fgamma4",strgamma4, -1,1);
1160 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE,1000);
1161 rf->SetParam(fgamma4,param->GetZWidth(), 1,0.2);
1162 rf->SetOffset(3*param->GetZSigma());
1165 TDirectory *savedir=gDirectory;
1166 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1168 AliFatal("Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !");
1171 prfinner->Read("prf_07504_Gati_056068_d02");
1172 //PH Set different names
1173 s=prfinner->GetGRF()->GetName();
1175 prfinner->GetGRF()->SetName(s.Data());
1177 prfouter1->Read("prf_10006_Gati_047051_d03");
1178 s=prfouter1->GetGRF()->GetName();
1180 prfouter1->GetGRF()->SetName(s.Data());
1182 prfouter2->Read("prf_15006_Gati_047051_d03");
1183 s=prfouter2->GetGRF()->GetName();
1185 prfouter2->GetGRF()->SetName(s.Data());
1190 param->SetInnerPRF(prfinner);
1191 param->SetOuter1PRF(prfouter1);
1192 param->SetOuter2PRF(prfouter2);
1193 param->SetTimeRF(rf);
1203 //__________________________________________________________________
1204 void AliTPC::Hits2Digits()
1207 // creates digits from hits
1209 if (!fTPCParam->IsGeoRead()){
1211 // read transformation matrices for gGeoManager
1213 fTPCParam->ReadGeoMatrices();
1216 fLoader->LoadHits("read");
1217 fLoader->LoadDigits("recreate");
1218 AliRunLoader* runLoader = fLoader->GetRunLoader();
1220 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1221 //PH runLoader->GetEvent(iEvent);
1222 Hits2Digits(iEvent);
1225 fLoader->UnloadHits();
1226 fLoader->UnloadDigits();
1228 //__________________________________________________________________
1229 void AliTPC::Hits2Digits(Int_t eventnumber)
1231 //----------------------------------------------------
1232 // Loop over all sectors for a single event
1233 //----------------------------------------------------
1234 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1235 rl->GetEvent(eventnumber);
1237 if (fLoader->TreeH() == 0x0) {
1238 if(fLoader->LoadHits()) {
1239 AliError("Can not load hits.");
1244 if (fLoader->TreeD() == 0x0 ) {
1245 fLoader->MakeTree("D");
1246 if (fLoader->TreeD() == 0x0 ) {
1247 AliError("Can not get TreeD");
1252 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1253 GenerNoise(500000); //create teble with noise
1255 //setup TPCDigitsArray
1257 if(GetDigitsArray()) delete GetDigitsArray();
1259 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1260 arr->SetClass("AliSimDigits");
1261 arr->Setup(fTPCParam);
1263 arr->MakeTree(fLoader->TreeD());
1264 SetDigitsArray(arr);
1266 fDigitsSwitch=0; // standard digits
1268 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1269 if (IsSectorActive(isec)) {
1270 AliDebug(1,Form("Hits2Digits","Sector %d is active.",isec));
1271 Hits2DigitsSector(isec);
1274 AliDebug(1,Form("Hits2Digits","Sector %d is NOT active.",isec));
1277 fLoader->WriteDigits("OVERWRITE");
1279 //this line prevents the crash in the similar one
1280 //on the beginning of this method
1281 //destructor attempts to reset the tree, which is deleted by the loader
1282 //need to be redesign
1283 if(GetDigitsArray()) delete GetDigitsArray();
1284 SetDigitsArray(0x0);
1288 //__________________________________________________________________
1289 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1292 //-----------------------------------------------------------
1293 // summable digits - 16 bit "ADC", no noise, no saturation
1294 //-----------------------------------------------------------
1296 //----------------------------------------------------
1297 // Loop over all sectors for a single event
1298 //----------------------------------------------------
1300 AliRunLoader* rl = fLoader->GetRunLoader();
1302 rl->GetEvent(eventnumber);
1303 if (fLoader->TreeH() == 0x0) {
1304 if(fLoader->LoadHits()) {
1305 AliError("Can not load hits.");
1312 if (fLoader->TreeS() == 0x0 ) {
1313 fLoader->MakeTree("S");
1316 if(fDefaults == 0) SetDefaults();
1318 GenerNoise(500000); //create table with noise
1319 //setup TPCDigitsArray
1321 if(GetDigitsArray()) delete GetDigitsArray();
1324 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1325 arr->SetClass("AliSimDigits");
1326 arr->Setup(fTPCParam);
1327 arr->MakeTree(fLoader->TreeS());
1329 SetDigitsArray(arr);
1331 fDigitsSwitch=1; // summable digits
1333 // set zero suppression to "0"
1335 fTPCParam->SetZeroSup(0);
1337 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1338 if (IsSectorActive(isec)) {
1339 Hits2DigitsSector(isec);
1342 fLoader->WriteSDigits("OVERWRITE");
1344 //this line prevents the crash in the similar one
1345 //on the beginning of this method
1346 //destructor attempts to reset the tree, which is deleted by the loader
1347 //need to be redesign
1348 if(GetDigitsArray()) delete GetDigitsArray();
1349 SetDigitsArray(0x0);
1351 //__________________________________________________________________
1353 void AliTPC::Hits2SDigits()
1356 //-----------------------------------------------------------
1357 // summable digits - 16 bit "ADC", no noise, no saturation
1358 //-----------------------------------------------------------
1359 if (0) fDebugStreamer = new TTreeSRedirector("TPCSimdebug.root");
1361 if (!fTPCParam->IsGeoRead()){
1363 // read transformation matrices for gGeoManager
1365 fTPCParam->ReadGeoMatrices();
1368 fLoader->LoadHits("read");
1369 fLoader->LoadSDigits("recreate");
1370 AliRunLoader* runLoader = fLoader->GetRunLoader();
1372 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1373 runLoader->GetEvent(iEvent);
1376 Hits2SDigits2(iEvent);
1379 fLoader->UnloadHits();
1380 fLoader->UnloadSDigits();
1381 if (fDebugStreamer) {
1382 delete fDebugStreamer;
1386 //_____________________________________________________________________________
1388 void AliTPC::Hits2DigitsSector(Int_t isec)
1390 //-------------------------------------------------------------------
1391 // TPC conversion from hits to digits.
1392 //-------------------------------------------------------------------
1394 //-----------------------------------------------------------------
1395 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1396 //-----------------------------------------------------------------
1398 //-------------------------------------------------------
1399 // Get the access to the track hits
1400 //-------------------------------------------------------
1402 // check if the parameters are set - important if one calls this method
1403 // directly, not from the Hits2Digits
1405 if(fDefaults == 0) SetDefaults();
1407 TTree *tH = fLoader->TreeH(); // pointer to the hits tree
1409 AliFatal("Can not find TreeH in folder");
1413 Stat_t ntracks = tH->GetEntries();
1419 Int_t nrows =fTPCParam->GetNRow(isec);
1421 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1423 MakeSector(isec,nrows,tH,ntracks,row);
1425 //--------------------------------------------------------
1426 // Digitize this sector, row by row
1427 // row[i] is the pointer to the TObjArray of TVectors,
1428 // each one containing electrons accepted on this
1429 // row, assigned into tracks
1430 //--------------------------------------------------------
1434 if (fDigitsArray->GetTree()==0) {
1435 AliFatal("Tree not set in fDigitsArray");
1438 for (i=0;i<nrows;i++){
1440 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1442 DigitizeRow(i,isec,row);
1444 fDigitsArray->StoreRow(isec,i);
1446 Int_t ndig = dig->GetDigitSize();
1449 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1452 fDigitsArray->ClearRow(isec,i);
1455 } // end of the sector digitization
1457 for(i=0;i<nrows+2;i++){
1462 delete [] row; // delete the array of pointers to TObjArray-s
1465 } // end of Hits2DigitsSector
1468 //_____________________________________________________________________________
1469 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1471 //-----------------------------------------------------------
1472 // Single row digitization, coupling from the neighbouring
1473 // rows taken into account
1474 //-----------------------------------------------------------
1476 //-----------------------------------------------------------------
1477 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1478 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1479 //-----------------------------------------------------------------
1481 Float_t zerosup = fTPCParam->GetZeroSup();
1482 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor();
1483 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
1484 AliTPCCalROC * gainROC = gainTPC->GetCalROC(isec); // pad gains per given sector
1485 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(isec); // noise per given sector
1488 fCurrentIndex[1]= isec;
1491 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1492 Int_t nofTbins = fTPCParam->GetMaxTBin();
1493 Int_t indexRange[4];
1495 // Integrated signal for this row
1496 // and a single track signal
1499 TMatrixF *m1 = new TMatrixF(0,nofPads,0,nofTbins); // integrated
1500 TMatrixF *m2 = new TMatrixF(0,nofPads,0,nofTbins); // single
1502 TMatrixF &total = *m1;
1504 // Array of pointers to the label-signal list
1506 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1507 Float_t **pList = new Float_t* [nofDigits];
1511 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1517 for (Int_t row= row1;row<=row2;row++){
1518 Int_t nTracks= rows[row]->GetEntries();
1519 for (i1=0;i1<nTracks;i1++){
1520 fCurrentIndex[2]= row;
1521 fCurrentIndex[3]=irow+1;
1523 m2->Zero(); // clear single track signal matrix
1524 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1525 GetList(trackLabel,nofPads,m2,indexRange,pList);
1527 else GetSignal(rows[row],i1,0,m1,indexRange);
1533 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1535 Float_t fzerosup = zerosup+0.5;
1536 for(Int_t it=0;it<nofTbins;it++){
1537 for(Int_t ip=0;ip<nofPads;ip++){
1539 Float_t q=total(ip,it);
1540 if(fDigitsSwitch == 0){
1541 Float_t gain = gainROC->GetValue(irow,ip); // get gain for given - pad-row pad
1542 Float_t noisePad = noiseROC->GetValue(irow,ip);
1545 q+=GetNoise()*noisePad;
1546 if(q <=fzerosup) continue; // do not fill zeros
1548 if(q >= fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat() - 1; // saturation
1553 if(q <= 0.) continue; // do not fill zeros
1554 if(q>2000.) q=2000.;
1560 // "real" signal or electronic noise (list = -1)?
1563 for(Int_t j1=0;j1<3;j1++){
1564 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1569 <A NAME="AliDigits"></A>
1570 using of AliDigits object
1573 dig->SetDigitFast((Short_t)q,it,ip);
1574 if (fDigitsArray->IsSimulated()) {
1575 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1576 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1577 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1580 } // end of loop over time buckets
1581 } // end of lop over pads
1587 // glitch filters if normal simulated digits
1589 if(!fDigitsSwitch) ((AliSimDigits*)dig)->GlitchFilter();
1591 // This row has been digitized, delete nonused stuff
1594 for(lp=0;lp<nofDigits;lp++){
1595 if(pList[lp]) delete [] pList[lp];
1603 } // end of DigitizeRow
1605 //_____________________________________________________________________________
1607 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1608 TMatrixF *m1, TMatrixF *m2,Int_t *indexRange)
1611 //---------------------------------------------------------------
1612 // Calculates 2-D signal (pad,time) for a single track,
1613 // returns a pointer to the signal matrix and the track label
1614 // No digitization is performed at this level!!!
1615 //---------------------------------------------------------------
1617 //-----------------------------------------------------------------
1618 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1619 // Modified: Marian Ivanov
1620 //-----------------------------------------------------------------
1624 tv = (TVector*)p1->At(ntr); // pointer to a track
1627 Float_t label = v(0);
1628 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1))/2;
1630 Int_t nElectrons = (tv->GetNrows()-1)/5;
1631 indexRange[0]=9999; // min pad
1632 indexRange[1]=-1; // max pad
1633 indexRange[2]=9999; //min time
1634 indexRange[3]=-1; // max time
1636 TMatrixF &signal = *m1;
1637 TMatrixF &total = *m2;
1639 // Loop over all electrons
1641 for(Int_t nel=0; nel<nElectrons; nel++){
1643 Float_t aval = v(idx+4);
1644 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1645 Float_t xyz[4]={v(idx+1),v(idx+2),v(idx+3),v(idx+5)};
1646 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
1648 Int_t *index = fTPCParam->GetResBin(0);
1649 Float_t *weight = & (fTPCParam->GetResWeight(0));
1651 if (n>0) for (Int_t i =0; i<n; i++){
1652 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1655 Int_t time=index[2];
1656 Float_t qweight = *(weight)*eltoadcfac;
1658 if (m1!=0) signal(pad,time)+=qweight;
1659 total(pad,time)+=qweight;
1660 if (indexRange[0]>pad) indexRange[0]=pad;
1661 if (indexRange[1]<pad) indexRange[1]=pad;
1662 if (indexRange[2]>time) indexRange[2]=time;
1663 if (indexRange[3]<time) indexRange[3]=time;
1670 } // end of loop over electrons
1672 return label; // returns track label when finished
1675 //_____________________________________________________________________________
1676 void AliTPC::GetList(Float_t label,Int_t np,TMatrixF *m,
1677 Int_t *indexRange, Float_t **pList)
1679 //----------------------------------------------------------------------
1680 // Updates the list of tracks contributing to digits for a given row
1681 //----------------------------------------------------------------------
1683 //-----------------------------------------------------------------
1684 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1685 //-----------------------------------------------------------------
1687 TMatrixF &signal = *m;
1689 // lop over nonzero digits
1691 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1692 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1695 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1697 if(signal(ip,it)<0.5) continue;
1699 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1701 if(!pList[globalIndex]){
1704 // Create new list (6 elements - 3 signals and 3 labels),
1707 pList[globalIndex] = new Float_t [6];
1711 *pList[globalIndex] = -1.;
1712 *(pList[globalIndex]+1) = -1.;
1713 *(pList[globalIndex]+2) = -1.;
1714 *(pList[globalIndex]+3) = -1.;
1715 *(pList[globalIndex]+4) = -1.;
1716 *(pList[globalIndex]+5) = -1.;
1718 *pList[globalIndex] = label;
1719 *(pList[globalIndex]+3) = signal(ip,it);
1723 // check the signal magnitude
1725 Float_t highest = *(pList[globalIndex]+3);
1726 Float_t middle = *(pList[globalIndex]+4);
1727 Float_t lowest = *(pList[globalIndex]+5);
1730 // compare the new signal with already existing list
1733 if(signal(ip,it)<lowest) continue; // neglect this track
1737 if (signal(ip,it)>highest){
1738 *(pList[globalIndex]+5) = middle;
1739 *(pList[globalIndex]+4) = highest;
1740 *(pList[globalIndex]+3) = signal(ip,it);
1742 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1743 *(pList[globalIndex]+1) = *pList[globalIndex];
1744 *pList[globalIndex] = label;
1746 else if (signal(ip,it)>middle){
1747 *(pList[globalIndex]+5) = middle;
1748 *(pList[globalIndex]+4) = signal(ip,it);
1750 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1751 *(pList[globalIndex]+1) = label;
1754 *(pList[globalIndex]+5) = signal(ip,it);
1755 *(pList[globalIndex]+2) = label;
1759 } // end of loop over pads
1760 } // end of loop over time bins
1763 //___________________________________________________________________
1764 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1765 Stat_t ntracks,TObjArray **row)
1768 //-----------------------------------------------------------------
1769 // Prepares the sector digitization, creates the vectors of
1770 // tracks for each row of this sector. The track vector
1771 // contains the track label and the position of electrons.
1772 //-----------------------------------------------------------------
1775 // The trasport of the electrons through TPC drift volume
1776 // Drift (drift velocity + velocity map(not yet implemented)))
1777 // Application of the random processes (diffusion, gas gain)
1778 // Systematic effects (ExB effect in drift volume + ROCs)
1781 // Loop over primary electrons:
1782 // Creation of the secondary electrons
1783 // Loop over electrons (primary+ secondaries)
1784 // Global coordinate frame:
1785 // 1. Skip electrons if attached
1786 // 2. ExB effect in drift volume
1787 // a.) Simulation calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1788 // b.) Reconstruction - calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1789 // 3. Generation of gas gain (Random - Exponential distribution)
1790 // 4. TransportElectron function (diffusion)
1792 // 5. Conversion to the local coordinate frame pad-row, pad, timebin
1793 // 6. Apply Time0 shift - AliTPCCalPad class
1794 // a.) Plus sign in simulation
1795 // b.) Minus sign in reconstruction
1798 //-----------------------------------------------------------------
1799 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1800 // Origin: Marian Ivanov, marian.ivanov@cern.ch
1801 //-----------------------------------------------------------------
1802 AliTPCcalibDB* const calib=AliTPCcalibDB::Instance();
1803 if (gAlice){ // Set correctly the magnetic field in the ExB calculation
1804 AliMagF * field = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField());
1806 calib->SetExBField(field->SolenoidField());
1810 Float_t gasgain = fTPCParam->GetGasGain();
1811 gasgain = gasgain/fGainFactor;
1815 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1818 if (fHitType>1) branch = TH->GetBranch("TPC2");
1819 else branch = TH->GetBranch("TPC");
1822 //----------------------------------------------
1823 // Create TObjArray-s, one for each row,
1824 // each TObjArray will store the TVectors
1825 // of electrons, one TVectors per each track.
1826 //----------------------------------------------
1828 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
1829 TVector **tracks = new TVector* [nrows+2]; //pointers to the track vectors
1831 for(i=0; i<nrows+2; i++){
1832 row[i] = new TObjArray;
1839 //--------------------------------------------------------------------
1840 // Loop over tracks, the "track" contains the full history
1841 //--------------------------------------------------------------------
1843 Int_t previousTrack,currentTrack;
1844 previousTrack = -1; // nothing to store so far!
1846 for(Int_t track=0;track<ntracks;track++){
1847 Bool_t isInSector=kTRUE;
1849 isInSector = TrackInVolume(isec,track);
1850 if (!isInSector) continue;
1852 branch->GetEntry(track); // get next track
1856 tpcHit = (AliTPChit*)FirstHit(-1);
1858 //--------------------------------------------------------------
1860 //--------------------------------------------------------------
1865 Int_t sector=tpcHit->fSector; // sector number
1867 tpcHit = (AliTPChit*) NextHit();
1871 // Remove hits which arrive before the TPC opening gate signal
1872 if(((fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
1873 /fTPCParam->GetDriftV()+tpcHit->Time())<fTPCParam->GetGateDelay()) {
1874 tpcHit = (AliTPChit*) NextHit();
1878 currentTrack = tpcHit->Track(); // track number
1880 if(currentTrack != previousTrack){
1882 // store already filled fTrack
1884 for(i=0;i<nrows+2;i++){
1885 if(previousTrack != -1){
1886 if(nofElectrons[i]>0){
1887 TVector &v = *tracks[i];
1888 v(0) = previousTrack;
1889 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
1890 row[i]->Add(tracks[i]);
1893 delete tracks[i]; // delete empty TVector
1899 tracks[i] = new TVector(601); // TVectors for the next fTrack
1901 } // end of loop over rows
1903 previousTrack=currentTrack; // update track label
1906 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1908 //---------------------------------------------------
1909 // Calculate the electron attachment probability
1910 //---------------------------------------------------
1913 Float_t time = 1.e6*(fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
1914 /fTPCParam->GetDriftV();
1916 Float_t attProb = fTPCParam->GetAttCoef()*
1917 fTPCParam->GetOxyCont()*time; // fraction!
1919 //-----------------------------------------------
1920 // Loop over electrons
1921 //-----------------------------------------------
1924 for(Int_t nel=0;nel<qI;nel++){
1925 // skip if electron lost due to the attachment
1926 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1931 Double_t dxyz0[3],dxyz1[3];
1932 dxyz0[0]=tpcHit->X();
1933 dxyz0[1]=tpcHit->Y();
1934 dxyz0[2]=tpcHit->Z();
1935 if (calib->GetExB()){
1936 calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1938 AliError("Not valid ExB calibration");
1939 dxyz1[0]=tpcHit->X();
1940 dxyz1[1]=tpcHit->Y();
1941 dxyz1[2]=tpcHit->Z();
1949 // protection for the nonphysical avalanche size (10**6 maximum)
1951 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
1952 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
1955 TransportElectron(xyz,index);
1957 Int_t padrow = fTPCParam->GetPadRow(xyz,index);
1959 // Add Time0 correction due unisochronity
1960 // xyz[0] - pad row coordinate
1961 // xyz[1] - pad coordinate
1962 // xyz[2] - is in now time bin coordinate system
1963 Float_t correction =0;
1964 if (calib->GetPadTime0()){
1965 if (!calib->GetPadTime0()->GetCalROC(isec)) continue;
1966 Int_t npads = fTPCParam->GetNPads(isec,padrow);
1967 // Int_t pad = TMath::Nint(xyz[1]+fTPCParam->GetNPads(isec,TMath::Nint(xyz[0]))*0.5);
1968 // pad numbering from -npads/2 .. npads/2-1
1969 Int_t pad = TMath::Nint(xyz[1]+npads/2);
1971 if (pad>=npads) pad=npads-1;
1972 correction = calib->GetPadTime0()->GetCalROC(isec)->GetValue(padrow,pad);
1973 // printf("%d\t%d\t%d\t%f\n",isec,padrow,pad,correction);
1974 if (fDebugStreamer){
1975 (*fDebugStreamer)<<"Time0"<<
1983 "cor="<<correction<<
1988 xyz[2]+=fTPCParam->GetNTBinsL1(); // adding Level 1 time bin offset
1990 // Electron track time (for pileup simulation)
1991 xyz[2]+=tpcHit->Time()/fTPCParam->GetTSample(); // adding time of flight
1995 // row 0 - cross talk from the innermost row
1996 // row fNRow+1 cross talk from the outermost row
1997 rowNumber = index[2]+1;
1998 //transform position to local digit coordinates
1999 //relative to nearest pad row
2000 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2002 if (isec <fTPCParam->GetNInnerSector()) {
2003 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2004 y1 = fTPCParam->GetYInner(rowNumber);
2007 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2008 y1 = fTPCParam->GetYOuter(rowNumber);
2010 // gain inefficiency at the wires edges - linear
2013 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.)); */
2015 nofElectrons[rowNumber]++;
2016 //----------------------------------
2017 // Expand vector if necessary
2018 //----------------------------------
2019 if(nofElectrons[rowNumber]>120){
2020 Int_t range = tracks[rowNumber]->GetNrows();
2021 if((nofElectrons[rowNumber])>(range-1)/5){
2023 tracks[rowNumber]->ResizeTo(range+500); // Add 100 electrons
2027 TVector &v = *tracks[rowNumber];
2028 Int_t idx = 5*nofElectrons[rowNumber]-4;
2029 Real_t * position = &(((TVector&)v)(idx)); //make code faster
2030 memcpy(position,xyz,5*sizeof(Float_t));
2032 } // end of loop over electrons
2034 tpcHit = (AliTPChit*)NextHit();
2036 } // end of loop over hits
2037 } // end of loop over tracks
2040 // store remaining track (the last one) if not empty
2043 for(i=0;i<nrows+2;i++){
2044 if(nofElectrons[i]>0){
2045 TVector &v = *tracks[i];
2046 v(0) = previousTrack;
2047 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
2048 row[i]->Add(tracks[i]);
2057 delete [] nofElectrons;
2059 } // end of MakeSector
2062 //_____________________________________________________________________________
2066 // Initialise TPC detector after definition of geometry
2068 AliDebug(1,"*********************************************");
2071 //_____________________________________________________________________________
2072 void AliTPC::ResetDigits()
2075 // Reset number of digits and the digits array for this detector
2078 if (fDigits) fDigits->Clear();
2083 //_____________________________________________________________________________
2084 void AliTPC::SetSens(Int_t sens)
2087 //-------------------------------------------------------------
2088 // Activates/deactivates the sensitive strips at the center of
2089 // the pad row -- this is for the space-point resolution calculations
2090 //-------------------------------------------------------------
2092 //-----------------------------------------------------------------
2093 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2094 //-----------------------------------------------------------------
2100 void AliTPC::SetSide(Float_t side=0.)
2102 // choice of the TPC side
2107 //_____________________________________________________________________________
2109 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2112 // electron transport taking into account:
2114 // 2.ExB at the wires
2115 // 3. nonisochronity
2117 // xyz and index must be already transformed to system 1
2120 fTPCParam->Transform1to2(xyz,index); // mis-alignment applied in this step
2123 Float_t driftl=xyz[2];
2124 if(driftl<0.01) driftl=0.01;
2125 driftl=TMath::Sqrt(driftl);
2126 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2127 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2128 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2129 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2130 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2134 if (fTPCParam->GetMWPCReadout()==kTRUE){
2135 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2136 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2138 //add nonisochronity (not implemented yet)
2144 //______________________________________________________________________
2145 AliTPChit::AliTPChit()
2157 //_____________________________________________________________________________
2158 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits)
2159 :AliHit(shunt,track),
2166 // Creates a TPC hit object
2177 //________________________________________________________________________
2178 // Additional code because of the AliTPCTrackHitsV2
2180 void AliTPC::MakeBranch(Option_t *option)
2183 // Create a new branch in the current Root Tree
2184 // The branch of fHits is automatically split
2185 // MI change 14.09.2000
2187 if (fHitType<2) return;
2188 char branchname[10];
2189 sprintf(branchname,"%s2",GetName());
2191 // Get the pointer to the header
2192 const char *cH = strstr(option,"H");
2194 if (fTrackHits && fLoader->TreeH() && cH && fHitType&4) {
2195 AliDebug(1,"Making branch for Type 4 Hits");
2196 fLoader->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2199 // if (fTrackHitsOld && fLoader->TreeH() && cH && fHitType&2) {
2200 // AliDebug(1,"Making branch for Type 2 Hits");
2201 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2202 // fLoader->TreeH(),fBufferSize,99);
2203 // fLoader->TreeH()->GetListOfBranches()->Add(branch);
2207 void AliTPC::SetTreeAddress()
2209 //Sets tree address for hits
2211 if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2212 AliDetector::SetTreeAddress();
2214 if (fHitType>1) SetTreeAddress2();
2217 void AliTPC::SetTreeAddress2()
2220 // Set branch address for the TrackHits Tree
2225 char branchname[20];
2226 sprintf(branchname,"%s2",GetName());
2228 // Branch address for hit tree
2229 TTree *treeH = fLoader->TreeH();
2230 if ((treeH)&&(fHitType&4)) {
2231 branch = treeH->GetBranch(branchname);
2233 branch->SetAddress(&fTrackHits);
2234 AliDebug(1,"fHitType&4 Setting");
2237 AliDebug(1,"fHitType&4 Failed (can not find branch)");
2240 // if ((treeH)&&(fHitType&2)) {
2241 // branch = treeH->GetBranch(branchname);
2243 // branch->SetAddress(&fTrackHitsOld);
2244 // AliDebug(1,"fHitType&2 Setting");
2247 // AliDebug(1,"fHitType&2 Failed (can not find branch)");
2251 void AliTPC::FinishPrimary()
2253 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2254 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2258 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2261 // add hit to the list
2265 int primary = gAlice->GetMCApp()->GetPrimary(track);
2266 gAlice->GetMCApp()->Particle(primary)->SetBit(kKeepBit);
2270 gAlice->GetMCApp()->FlagTrack(track);
2272 if (fTrackHits && fHitType&4)
2273 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2274 hits[1],hits[2],(Int_t)hits[3],hits[4]);
2275 // if (fTrackHitsOld &&fHitType&2 )
2276 // fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2277 // hits[1],hits[2],(Int_t)hits[3]);
2281 void AliTPC::ResetHits()
2283 if (fHitType&1) AliDetector::ResetHits();
2284 if (fHitType>1) ResetHits2();
2287 void AliTPC::ResetHits2()
2291 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2292 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2296 AliHit* AliTPC::FirstHit(Int_t track)
2298 if (fHitType>1) return FirstHit2(track);
2299 return AliDetector::FirstHit(track);
2301 AliHit* AliTPC::NextHit()
2306 if (fHitType>1) return NextHit2();
2308 return AliDetector::NextHit();
2311 AliHit* AliTPC::FirstHit2(Int_t track)
2314 // Initialise the hit iterator
2315 // Return the address of the first hit for track
2316 // If track>=0 the track is read from disk
2317 // while if track<0 the first hit of the current
2318 // track is returned
2321 gAlice->GetMCApp()->ResetHits();
2322 fLoader->TreeH()->GetEvent(track);
2325 if (fTrackHits && fHitType&4) {
2326 fTrackHits->First();
2327 return fTrackHits->GetHit();
2329 // if (fTrackHitsOld && fHitType&2) {
2330 // fTrackHitsOld->First();
2331 // return fTrackHitsOld->GetHit();
2337 AliHit* AliTPC::NextHit2()
2340 //Return the next hit for the current track
2343 // if (fTrackHitsOld && fHitType&2) {
2344 // fTrackHitsOld->Next();
2345 // return fTrackHitsOld->GetHit();
2349 return fTrackHits->GetHit();
2355 void AliTPC::RemapTrackHitIDs(Int_t *map)
2360 if (!fTrackHits) return;
2362 // if (fTrackHitsOld && fHitType&2){
2363 // AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2364 // for (UInt_t i=0;i<arr->GetSize();i++){
2365 // AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2366 // info->fTrackID = map[info->fTrackID];
2369 // if (fTrackHitsOld && fHitType&4){
2370 if (fTrackHits && fHitType&4){
2371 TClonesArray * arr = fTrackHits->GetArray();;
2372 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2373 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2374 info->SetTrackID(map[info->GetTrackID()]);
2379 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2381 //return bool information - is track in given volume
2382 //load only part of the track information
2383 //return true if current track is in volume
2386 // if (fTrackHitsOld && fHitType&2) {
2387 // TBranch * br = fLoader->TreeH()->GetBranch("fTrackHitsInfo");
2388 // br->GetEvent(track);
2389 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2390 // for (UInt_t j=0;j<ar->GetSize();j++){
2391 // if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2395 if (fTrackHits && fHitType&4) {
2396 TBranch * br1 = fLoader->TreeH()->GetBranch("fVolumes");
2397 TBranch * br2 = fLoader->TreeH()->GetBranch("fNVolumes");
2398 br2->GetEvent(track);
2399 br1->GetEvent(track);
2400 Int_t *volumes = fTrackHits->GetVolumes();
2401 Int_t nvolumes = fTrackHits->GetNVolumes();
2402 if (!volumes && nvolumes>0) {
2403 AliWarning(Form("Problematic track\t%d\t%d",track,nvolumes));
2406 for (Int_t j=0;j<nvolumes; j++)
2407 if (volumes[j]==id) return kTRUE;
2411 TBranch * br = fLoader->TreeH()->GetBranch("fSector");
2412 br->GetEvent(track);
2413 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2414 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2422 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
2425 fLoader = new AliTPCLoader(GetName(),topfoldername);
2429 ////////////////////////////////////////////////////////////////////////
2430 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
2432 // load TPC paarmeters from a given file or create new if the object
2433 // is not found there
2434 // 12/05/2003 This method should be moved to the AliTPCLoader
2435 // and one has to decide where to store the TPC parameters
2438 sprintf(paramName,"75x40_100x60_150x60");
2439 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
2441 AliDebugClass(1,Form("TPC parameters %s found.",paramName));
2443 AliWarningClass("TPC parameters not found. Create new (they may be incorrect)");
2444 //paramTPC = new AliTPCParamSR;
2445 paramTPC = AliTPCcalibDB::Instance()->GetParameters();
2446 if (!paramTPC->IsGeoRead()){
2448 // read transformation matrices for gGeoManager
2450 paramTPC->ReadGeoMatrices();
2456 // the older version of parameters can be accessed with this code.
2457 // In some cases, we have old parameters saved in the file but
2458 // digits were created with new parameters, it can be distinguish
2459 // by the name of TPC TreeD. The code here is just for the case
2460 // we would need to compare with old data, uncomment it if needed.
2462 // char paramName[50];
2463 // sprintf(paramName,"75x40_100x60");
2464 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
2466 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2468 // sprintf(paramName,"75x40_100x60_150x60");
2469 // paramTPC=(AliTPCParam*)in->Get(paramName);
2471 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2473 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
2475 // paramTPC = new AliTPCParamSR;