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),
109 // Default constructor
112 for(Int_t i=0;i<4;i++) fCurrentIndex[i]=0;
114 // fTrackHitsOld = 0;
115 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
116 fHitType = 4; // ROOT containers
118 fHitType = 2; //default CONTAINERS - based on ROOT structure
122 //_____________________________________________________________________________
123 AliTPC::AliTPC(const char *name, const char *title)
124 : AliDetector(name,title),
134 fPrimaryIonisation(0),
146 // Standard constructor
150 // Initialise arrays of hits and digits
151 fHits = new TClonesArray("AliTPChit", 176);
152 gAlice->GetMCApp()->AddHitList(fHits);
154 fTrackHits = new AliTPCTrackHitsV2;
155 fTrackHits->SetHitPrecision(0.002);
156 fTrackHits->SetStepPrecision(0.003);
157 fTrackHits->SetMaxDistance(100);
159 //fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
160 //fTrackHitsOld->SetHitPrecision(0.002);
161 //fTrackHitsOld->SetStepPrecision(0.003);
162 //fTrackHitsOld->SetMaxDistance(100);
165 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
166 fHitType = 4; // ROOT containers
171 for(Int_t i=0;i<4;i++) fCurrentIndex[i]=0;
176 // Initialise color attributes
177 //PH SetMarkerColor(kYellow);
180 // Set TPC parameters
184 if (!strcmp(title,"Ne-CO2") || !strcmp(title,"Ne-CO2-N") || !strcmp(title,"Default") ) {
185 fTPCParam = AliTPCcalibDB::Instance()->GetParameters();
188 AliWarning("In Config.C you must set non-default parameters.");
193 void AliTPC::CreateDebugStremer(){
195 // Create Debug streamer to check simulation
197 fDebugStreamer = new TTreeSRedirector("TPCSimdebug.root");
199 //_____________________________________________________________________________
210 delete fTrackHits; //MI 15.09.2000
211 // delete fTrackHitsOld; //MI 10.12.2001
214 delete [] fNoiseTable;
215 delete [] fActiveSectors;
216 if (fDebugStreamer) delete fDebugStreamer;
219 //_____________________________________________________________________________
220 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
223 // Add a hit to the list
226 TClonesArray &lhits = *fHits;
227 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
230 AddHit2(track,vol,hits);
233 //_____________________________________________________________________________
234 void AliTPC::CreateMaterials()
236 //-----------------------------------------------
237 // Create Materials for for TPC simulations
238 //-----------------------------------------------
240 //-----------------------------------------------------------------
241 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
242 //-----------------------------------------------------------------
244 Int_t iSXFLD=((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ();
245 Float_t sXMGMX=((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();
247 Float_t amat[5]; // atomic numbers
248 Float_t zmat[5]; // z
249 Float_t wmat[5]; // proportions
255 //***************** Gases *************************
258 //--------------------------------------------------------------
259 // gases - air and CO2
260 //--------------------------------------------------------------
276 AliMixture(10,"CO2",amat,zmat,density,2,wmat);
291 AliMixture(11,"Air",amat,zmat,density,2,wmat);
293 //----------------------------------------------------------------
295 //----------------------------------------------------------------
298 // Drift gases 1 - nonsensitive, 2 - sensitive, 3 - for Kr
299 // Composition by % of volume, values at 20deg and 1 atm.
301 // get the geometry title - defined in Config.C
304 TString title(GetTitle());
316 if(title == TString("Ne-CO2") || title == TString("Default")){
322 AliMixture(12,"Ne-CO2-1",amat,zmat,density,3,wmat);
323 AliMixture(13,"Ne-CO2-2",amat,zmat,density,3,wmat);
324 AliMixture(35,"Ne-CO2-3",amat,zmat,density,3,wmat);
326 else if (title == TString("Ne-CO2-N")){
335 AliMixture(12,"Ne-CO2-N-1",amat,zmat,density,4,wmat);
336 AliMixture(13,"Ne-CO2-N-2",amat,zmat,density,4,wmat);
337 AliMixture(30,"Ne-CO2-N-3",amat,zmat,density,4,wmat);
343 //----------------------------------------------------------------------
345 //----------------------------------------------------------------------
367 AliMixture(14,"Kevlar",amat,zmat,density,-4,wmat);
388 AliMixture(15,"NOMEX",amat,zmat,density,-4,wmat);
406 AliMixture(16,"Makrolon",amat,zmat,density,-3,wmat);
424 AliMixture(17, "Tedlar",amat,zmat,density,-3,wmat);
442 AliMixture(18, "Mylar",amat,zmat,density,-3,wmat);
443 // material for "prepregs"
444 // Epoxy - C14 H20 O3
447 // prepreg1 60% C-fiber, 40% epoxy (vol)
462 AliMixture(19, "Prepreg1",amat,zmat,density,3,wmat);
464 //prepreg2 60% glass-fiber, 40% epoxy
483 AliMixture(20, "Prepreg2",amat,zmat,density,4,wmat);
485 //prepreg3 50% glass-fiber, 50% epoxy
504 AliMixture(21, "Prepreg3",amat,zmat,density,4,wmat);
506 // G10 60% SiO2 40% epoxy
525 AliMixture(22, "G10",amat,zmat,density,4,wmat);
534 AliMaterial(23,"Al",amat[0],zmat[0],density,999.,999.);
536 // Si (for electronics
543 AliMaterial(24,"Si",amat[0],zmat[0],density,999.,999.);
552 AliMaterial(25,"Cu",amat[0],zmat[0],density,999.,999.);
570 AliMixture(33,"Brass",amat,zmat,density,2,wmat);
572 // Epoxy - C14 H20 O3
588 AliMixture(26,"Epoxy",amat,zmat,density,-3,wmat);
590 // epoxy film - 90% epoxy, 10% glass fiber
610 AliMixture(34, "Epoxy-film",amat,zmat,density,4,wmat);
628 AliMixture(27,"Plexiglas",amat,zmat,density,-3,wmat);
636 AliMaterial(28,"C",amat[0],zmat[0],density,999.,999.);
638 // Fe (steel for the inner heat screen)
646 AliMaterial(29,"Fe",amat[0],zmat[0],density,999.,999.);
648 // Peek - (C6H4-O-OC6H4-O-C6H4-CO)n
663 AliMixture(30,"Peek",amat,zmat,density,-3,wmat);
678 AliMixture(31,"Alumina",amat,zmat,density,-2,wmat);
697 AliMixture(32,"Water",amat,zmat,density,-2,wmat);
700 //----------------------------------------------------------
701 // tracking media for gases
702 //----------------------------------------------------------
704 AliMedium(0, "Air", 11, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
705 AliMedium(1, "DriftGas1", 12, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
706 AliMedium(2, "DriftGas2", 13, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
707 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
708 AliMedium(20, "DriftGas3", 35, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
709 //-----------------------------------------------------------
710 // tracking media for solids
711 //-----------------------------------------------------------
713 AliMedium(4,"Al",23,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
714 AliMedium(5,"Kevlar",14,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
715 AliMedium(6,"Nomex",15,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
716 AliMedium(7,"Makrolon",16,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
717 AliMedium(8,"Mylar",18,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
718 AliMedium(9,"Tedlar",17,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
720 AliMedium(10,"Prepreg1",19,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
721 AliMedium(11,"Prepreg2",20,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
722 AliMedium(12,"Prepreg3",21,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
723 AliMedium(13,"Epoxy",26,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
725 AliMedium(14,"Cu",25,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
726 AliMedium(15,"Si",24,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
727 AliMedium(16,"G10",22,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
728 AliMedium(17,"Plexiglas",27,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
729 AliMedium(18,"Steel",29,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
730 AliMedium(19,"Peek",30,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
731 AliMedium(21,"Alumina",31,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
732 AliMedium(22,"Water",32,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
733 AliMedium(23,"Brass",33,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
734 AliMedium(24,"Epoxyfm",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
737 void AliTPC::GenerNoise(Int_t tablesize)
740 //Generate table with noise
747 if (fNoiseTable) delete[] fNoiseTable;
748 fNoiseTable = new Float_t[tablesize];
749 fNoiseDepth = tablesize;
750 fCurrentNoise =0; //!index of the noise in the noise table
752 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
753 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
756 Float_t AliTPC::GetNoise()
758 // get noise from table
759 // if ((fCurrentNoise%10)==0)
760 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
761 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
762 return fNoiseTable[fCurrentNoise++];
763 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
767 Bool_t AliTPC::IsSectorActive(Int_t sec) const
770 // check if the sector is active
771 if (!fActiveSectors) return kTRUE;
772 else return fActiveSectors[sec];
775 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
777 // activate interesting sectors
778 SetTreeAddress();//just for security
779 if (!fActiveSectors) fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
780 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
781 for (Int_t i=0;i<n;i++)
782 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
786 void AliTPC::SetActiveSectors(Int_t flag)
789 // activate sectors which were hitted by tracks
791 SetTreeAddress();//just for security
792 if (fHitType==0) return; // if Clones hit - not short volume ID information
793 if (!fActiveSectors) fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
795 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
798 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
799 //TBranch * branch=0;
800 if (fLoader->TreeH() == 0x0)
802 AliFatal("Can not find TreeH in folder");
805 //if (fHitType>1) branch = fLoader->TreeH()->GetBranch("TPC2");
806 if (fHitType>1) fLoader->TreeH()->GetBranch("TPC2");
807 //else branch = fLoader->TreeH()->GetBranch("TPC");
808 else fLoader->TreeH()->GetBranch("TPC");
809 Stat_t ntracks = fLoader->TreeH()->GetEntries();
810 // loop over all hits
811 AliDebug(1,Form("Got %d tracks", (Int_t) ntracks));
813 for(Int_t track=0;track<ntracks;track++) {
816 if (fTrackHits && fHitType&4) {
817 TBranch * br1 = fLoader->TreeH()->GetBranch("fVolumes");
818 TBranch * br2 = fLoader->TreeH()->GetBranch("fNVolumes");
819 br1->GetEvent(track);
820 br2->GetEvent(track);
821 Int_t *volumes = fTrackHits->GetVolumes();
822 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++) {
823 if (volumes[j]>-1 && volumes[j]<fTPCParam->GetNSector()) {
824 fActiveSectors[volumes[j]]=kTRUE;
827 AliError(Form("Volume %d -> sector number %d is outside (0..%d)",
830 fTPCParam->GetNSector()));
836 // if (fTrackHitsOld && fHitType&2) {
837 // TBranch * br = fLoader->TreeH()->GetBranch("fTrackHitsInfo");
838 // br->GetEvent(track);
839 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
840 // for (UInt_t j=0;j<ar->GetSize();j++){
841 // fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
850 //_____________________________________________________________________________
851 void AliTPC::Digits2Raw()
853 // convert digits of the current event to raw data
855 static const Int_t kThreshold = 0;
857 fLoader->LoadDigits();
858 TTree* digits = fLoader->TreeD();
860 AliError("No digits tree");
866 AliSimDigits* digrow = &digarr;
867 digits->GetBranch("Segment")->SetAddress(&digrow);
869 const char* fileName = "AliTPCDDL.dat";
870 AliTPCBuffer* buffer = new AliTPCBuffer(fileName);
874 // 2: txt files with digits
875 //BE CAREFUL, verbose level 2 MUST be used only for debugging and
876 //it is highly suggested to use this mode only for debugging digits files
877 //reasonably small, because otherwise the size of the txt files can reach
878 //quickly several MB wasting time and disk space.
879 buffer->SetVerbose(0);
881 Int_t nEntries = Int_t(digits->GetEntries());
882 Int_t previousSector = -1;
884 for (Int_t i = 0; i < nEntries; i++) {
887 fTPCParam->AdjustSectorRow(digarr.GetID(), sector, row);
888 if(previousSector != sector) {
890 previousSector = sector;
893 if (sector < 36) { //inner sector [0;35]
895 //the whole row is written into the output file
896 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
897 sector, subSector, row);
899 //only the pads in the range [37;48] are written into the output file
900 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 1,
901 sector, subSector, row);
903 //only the pads outside the range [37;48] are written into the output file
904 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 2,
905 sector, subSector, row);
908 } else { //outer sector [36;71]
909 if (row == 54) subSector = 2;
910 if ((row != 27) && (row != 76)) {
911 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
912 sector, subSector, row);
913 } else if (row == 27) {
914 //only the pads outside the range [43;46] are written into the output file
915 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 2,
916 sector, subSector, row);
918 //only the pads in the range [43;46] are written into the output file
919 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 1,
920 sector, subSector, row);
921 } else if (row == 76) {
922 //only the pads outside the range [33;88] are written into the output file
923 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 2,
924 sector, subSector, row);
926 //only the pads in the range [33;88] are written into the output file
927 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 1,
928 sector, subSector, row);
934 fLoader->UnloadDigits();
936 AliTPCDDLRawData rawWriter;
937 rawWriter.SetVerbose(0);
939 rawWriter.RawData(fileName);
940 gSystem->Unlink(fileName);
945 //_____________________________________________________________________________
946 Bool_t AliTPC::Raw2SDigits(AliRawReader* rawReader){
947 // Converts the TPC raw data into summable digits
948 // The method is used for merging simulated and
950 if (fLoader->TreeS() == 0x0 ) {
951 fLoader->MakeTree("S");
954 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
956 //setup TPCDigitsArray
957 if(GetDigitsArray()) delete GetDigitsArray();
959 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
960 arr->SetClass("AliSimDigits");
961 arr->Setup(fTPCParam);
962 arr->MakeTree(fLoader->TreeS());
966 // set zero suppression to "0"
967 fTPCParam->SetZeroSup(0);
970 const Int_t kmaxTime = fTPCParam->GetMaxTBin();
971 const Int_t kNIS = fTPCParam->GetNInnerSector();
972 const Int_t kNOS = fTPCParam->GetNOuterSector();
973 const Int_t kNS = kNIS + kNOS;
976 AliTPCROC * roc = AliTPCROC::Instance();
977 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
978 Int_t nPadsMax = roc->GetNPads(roc->GetNSector()-1,nRowsMax-1);
979 Short_t** allBins = new Short_t*[nRowsMax];
980 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
981 Int_t maxBin = kmaxTime*nPadsMax;
982 allBins[iRow] = new Short_t[maxBin];
983 memset(allBins[iRow],0,sizeof(Short_t)*maxBin);
986 for(Int_t iSector = 0; iSector < kNS; iSector++) {
988 Int_t nRows = fTPCParam->GetNRow(iSector);
989 Int_t nDDLs = 0, indexDDL = 0;
990 if (iSector < kNIS) {
992 indexDDL = iSector * 2;
996 indexDDL = (iSector-kNIS) * 4 + kNIS * 2;
999 // Load the raw data for corresponding DDLs
1002 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();
1003 AliTPCRawStreamV3 input(rawReader,(AliAltroMapping**)mapping);
1004 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
1007 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
1008 Int_t maxBin = kmaxTime*nPadsMax;
1009 memset(allBins[iRow],0,sizeof(Short_t)*maxBin);
1012 // Begin loop over altro data
1013 while (input.NextDDL()) {
1015 if (input.GetSector() != iSector)
1016 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",iSector,input.GetSector()));
1019 while ( input.NextChannel() ) {
1021 Int_t iRow = input.GetRow();
1022 if (iRow < 0 || iRow >= nRows)
1023 AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
1024 iRow, 0, nRows -1));
1025 Int_t iPad = input.GetPad();
1027 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
1029 if (iPad < 0 || iPad >= maxPad)
1030 AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
1031 iPad, 0, maxPad -1));
1034 while ( input.NextBunch() ){
1035 Int_t startTbin = (Int_t)input.GetStartTimeBin();
1036 Int_t bunchlength = (Int_t)input.GetBunchLength();
1037 const UShort_t *sig = input.GetSignals();
1038 for (Int_t iTime = 0; iTime<bunchlength; iTime++){
1039 Int_t iTimeBin=startTbin-iTime;
1040 if ( iTimeBin < 0 || iTimeBin >= kmaxTime) {
1042 //AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
1043 // iTimeBin, 0, kmaxTime -1));
1046 Int_t maxBin = kmaxTime*maxPad;
1047 if (((iPad*kmaxTime+iTimeBin) >= maxBin) ||
1048 ((iPad*kmaxTime+iTimeBin) < 0))
1049 AliFatal(Form("Index outside the allowed range"
1050 " Sector=%d Row=%d Pad=%d Timebin=%d"
1051 " (Max.index=%d)",iSector,iRow,iPad,iTimeBin,maxBin));
1052 allBins[iRow][iPad*kmaxTime+iTimeBin] = sig[iTime];
1055 } // End loop over altro data
1058 // Now fill the digits array
1059 if (fDigitsArray->GetTree()==0) {
1060 AliFatal("Tree not set in fDigitsArray");
1063 for (Int_t iRow = 0; iRow < nRows; iRow++) {
1064 AliDigits * dig = fDigitsArray->CreateRow(iSector,iRow);
1065 Int_t maxPad = fTPCParam->GetNPads(iSector,iRow);
1066 for(Int_t iPad = 0; iPad < maxPad; iPad++) {
1067 for(Int_t iTimeBin = 0; iTimeBin < kmaxTime; iTimeBin++) {
1068 Short_t q = allBins[iRow][iPad*kmaxTime + iTimeBin];
1069 if (q <= 0) continue;
1071 dig->SetDigitFast((Short_t)q,iTimeBin,iPad);
1072 ((AliSimDigits*)dig)->SetTrackIDFast( 3141593, iTimeBin,iPad,0);
1075 fDigitsArray->StoreRow(iSector,iRow);
1076 Int_t ndig = dig->GetDigitSize();
1079 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1080 iSector,iRow,ndig));
1082 fDigitsArray->ClearRow(iSector,iRow);
1084 } // end of the sector digitization
1086 // get LHC clock phase from the digits tree
1088 TParameter<float> *ph;
1090 TTree *digtree = fLoader->TreeD();
1092 if(digtree){ // if TreeD exists
1093 ph = (TParameter<float>*)digtree->GetUserInfo()->FindObject("lhcphase0");
1094 phase = ph->GetVal();
1096 else{ //TreeD does not exist
1100 // store lhc clock phase in S-digits tree
1102 fLoader->TreeS()->GetUserInfo()->Add(new TParameter<float>("lhcphase0",phase));
1104 fLoader->WriteSDigits("OVERWRITE");
1106 if(GetDigitsArray()) delete GetDigitsArray();
1107 SetDigitsArray(0x0);
1110 for (Int_t iRow = 0; iRow < nRowsMax; iRow++)
1111 delete [] allBins[iRow];
1117 //______________________________________________________________________
1118 AliDigitizer* AliTPC::CreateDigitizer(AliDigitizationInput* digInput) const
1120 return new AliTPCDigitizer(digInput);
1123 void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
1125 //create digits from summable digits
1126 GenerNoise(500000); //create teble with noise
1128 //conect tree with sSDigits
1129 TTree *t = fLoader->TreeS();
1132 fLoader->LoadSDigits("READ");
1133 t = fLoader->TreeS();
1135 AliError("Can not get input TreeS");
1140 if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
1142 AliSimDigits digarr, *dummy=&digarr;
1143 TBranch* sdb = t->GetBranch("Segment");
1145 AliError("Can not find branch with segments in TreeS.");
1149 sdb->SetAddress(&dummy);
1151 Stat_t nentries = t->GetEntries();
1153 // set zero suppression
1155 fTPCParam->SetZeroSup(2);
1157 // get zero suppression
1159 Int_t zerosup = fTPCParam->GetZeroSup();
1161 //make tree with digits
1163 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1164 arr->SetClass("AliSimDigits");
1165 arr->Setup(fTPCParam);
1166 arr->MakeTree(fLoader->TreeD());
1168 AliTPCParam * par = fTPCParam;
1170 //Loop over segments of the TPC
1172 for (Int_t n=0; n<nentries; n++) {
1175 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1176 AliWarning(Form("Invalid segment ID ! %d",digarr.GetID()));
1179 if (!IsSectorActive(sec)) continue;
1181 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1182 Int_t nrows = digrow->GetNRows();
1183 Int_t ncols = digrow->GetNCols();
1185 digrow->ExpandBuffer();
1186 digarr.ExpandBuffer();
1187 digrow->ExpandTrackBuffer();
1188 digarr.ExpandTrackBuffer();
1191 Short_t * pamp0 = digarr.GetDigits();
1192 Int_t * ptracks0 = digarr.GetTracks();
1193 Short_t * pamp1 = digrow->GetDigits();
1194 Int_t * ptracks1 = digrow->GetTracks();
1195 Int_t nelems =nrows*ncols;
1196 Int_t saturation = fTPCParam->GetADCSat() - 1;
1197 //use internal structure of the AliDigits - for speed reason
1198 //if you cahnge implementation
1199 //of the Alidigits - it must be rewriten -
1200 for (Int_t i= 0; i<nelems; i++){
1201 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1203 if (q>saturation) q=saturation;
1206 ptracks1[0]=ptracks0[0];
1207 ptracks1[nelems]=ptracks0[nelems];
1208 ptracks1[2*nelems]=ptracks0[2*nelems];
1216 arr->StoreRow(sec,row);
1217 arr->ClearRow(sec,row);
1222 fLoader->WriteDigits("OVERWRITE");
1226 //__________________________________________________________________
1227 void AliTPC::SetDefaults(){
1229 // setting the defaults
1232 // Set response functions
1235 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1237 //AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1238 //gDirectory->Get("75x40_100x60");
1239 AliTPCParamSR *param = (AliTPCParamSR*)AliTPCcalibDB::Instance()->GetParameters();
1241 AliFatal("No TPC parameters found");
1244 if (!param->IsGeoRead()){
1246 // read transformation matrices for gGeoManager
1248 param->ReadGeoMatrices();
1253 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1254 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1255 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1258 //AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1259 //rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1260 //rf->SetOffset(3*param->GetZSigma());
1265 char strgamma4[1000];
1266 //sprintf(strgamma4,"AliTPCRF1D::Gamma4((x-0.135+%f)*%f,55,160)",3*param->GetZSigma(), 1000000000*param->GetTSample()/param->GetZWidth());
1268 snprintf(strgamma4,1000,"AliTPCRF1D::Gamma4((x-0.135+%f)*%f,55,160)",3*param->GetZSigma(), 1000000000*param->GetTSample()/param->GetZWidth());
1269 TF1 * fgamma4 = new TF1("fgamma4",strgamma4, -1,1);
1270 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE,1000);
1271 rf->SetParam(fgamma4,param->GetZWidth(), 1,0.2);
1272 rf->SetOffset(3*param->GetZSigma());
1274 TDirectory *savedir=gDirectory;
1277 printf ("TPC MWPC readout\n");
1278 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1280 AliFatal("Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !");
1283 prfinner->Read("prf_07504_Gati_056068_d02");
1284 //PH Set different names
1285 s=prfinner->GetGRF()->GetName();
1287 prfinner->GetGRF()->SetName(s.Data());
1289 prfouter1->Read("prf_10006_Gati_047051_d03");
1290 s=prfouter1->GetGRF()->GetName();
1292 prfouter1->GetGRF()->SetName(s.Data());
1294 prfouter2->Read("prf_15006_Gati_047051_d03");
1295 s=prfouter2->GetGRF()->GetName();
1297 prfouter2->GetGRF()->SetName(s.Data());
1302 printf ("TPC GEM readout\n");
1303 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2dGEM.root");
1305 AliFatal("Can't open $ALICE_ROOT/TPC/AliTPCprf2dGEM.root !");
1308 prfinner->Read("prf0");
1309 //PH Set different names
1310 s=prfinner->GetGRF()->GetName();
1312 prfinner->GetGRF()->SetName(s.Data());
1314 prfouter1->Read("prf1");
1315 s=prfouter1->GetGRF()->GetName();
1317 prfouter1->GetGRF()->SetName(s.Data());
1319 prfouter2->Read("prf2");
1320 s=prfouter2->GetGRF()->GetName();
1322 prfouter2->GetGRF()->SetName(s.Data());
1327 param->SetInnerPRF(prfinner);
1328 param->SetOuter1PRF(prfouter1);
1329 param->SetOuter2PRF(prfouter2);
1330 param->SetTimeRF(rf);
1340 //__________________________________________________________________
1341 void AliTPC::Hits2Digits()
1344 // creates digits from hits
1346 if (!fTPCParam->IsGeoRead()){
1348 // read transformation matrices for gGeoManager
1350 fTPCParam->ReadGeoMatrices();
1353 fLoader->LoadHits("read");
1354 fLoader->LoadDigits("recreate");
1355 AliRunLoader* runLoader = fLoader->GetRunLoader();
1357 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1358 //PH runLoader->GetEvent(iEvent);
1359 Hits2Digits(iEvent);
1362 fLoader->UnloadHits();
1363 fLoader->UnloadDigits();
1365 //__________________________________________________________________
1366 void AliTPC::Hits2Digits(Int_t eventnumber)
1368 //----------------------------------------------------
1369 // Loop over all sectors for a single event
1370 //----------------------------------------------------
1371 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1372 rl->GetEvent(eventnumber);
1374 if (fLoader->TreeH() == 0x0) {
1375 if(fLoader->LoadHits()) {
1376 AliError("Can not load hits.");
1381 if (fLoader->TreeD() == 0x0 ) {
1382 fLoader->MakeTree("D");
1383 if (fLoader->TreeD() == 0x0 ) {
1384 AliError("Can not get TreeD");
1389 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1390 GenerNoise(500000); //create teble with noise
1392 //setup TPCDigitsArray
1394 if(GetDigitsArray()) delete GetDigitsArray();
1396 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1397 arr->SetClass("AliSimDigits");
1398 arr->Setup(fTPCParam);
1400 arr->MakeTree(fLoader->TreeD());
1401 SetDigitsArray(arr);
1403 fDigitsSwitch=0; // standard digits
1404 // here LHC clock phase
1406 switch (fLHCclockPhaseSw){
1413 lhcph = (Int_t)(gRandom->Rndm()/0.25);
1417 // not implemented yet
1420 // adding phase to the TreeD user info
1421 fLoader->TreeD()->GetUserInfo()->Add(new TParameter<float>("lhcphase0",lhcph));
1423 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1424 if (IsSectorActive(isec)) {
1425 AliDebug(1,Form("Hits2Digits: Sector %d is active.",isec));
1426 Hits2DigitsSector(isec);
1429 AliDebug(1,Form("Hits2Digits: Sector %d is NOT active.",isec));
1432 fLoader->WriteDigits("OVERWRITE");
1434 //this line prevents the crash in the similar one
1435 //on the beginning of this method
1436 //destructor attempts to reset the tree, which is deleted by the loader
1437 //need to be redesign
1438 if(GetDigitsArray()) delete GetDigitsArray();
1439 SetDigitsArray(0x0);
1443 //__________________________________________________________________
1444 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1447 //-----------------------------------------------------------
1448 // summable digits - 16 bit "ADC", no noise, no saturation
1449 //-----------------------------------------------------------
1451 //----------------------------------------------------
1452 // Loop over all sectors for a single event
1453 //----------------------------------------------------
1455 AliRunLoader* rl = fLoader->GetRunLoader();
1457 rl->GetEvent(eventnumber);
1458 if (fLoader->TreeH() == 0x0) {
1459 if(fLoader->LoadHits()) {
1460 AliError("Can not load hits.");
1467 if (fLoader->TreeS() == 0x0 ) {
1468 fLoader->MakeTree("S");
1471 if(fDefaults == 0) SetDefaults();
1473 GenerNoise(500000); //create table with noise
1474 //setup TPCDigitsArray
1476 if(GetDigitsArray()) delete GetDigitsArray();
1479 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1480 arr->SetClass("AliSimDigits");
1481 arr->Setup(fTPCParam);
1482 arr->MakeTree(fLoader->TreeS());
1484 SetDigitsArray(arr);
1486 fDigitsSwitch=1; // summable digits
1488 // set zero suppression to "0"
1489 // here LHC clock phase
1491 switch (fLHCclockPhaseSw){
1498 lhcph = (Int_t)(gRandom->Rndm()/0.25);
1502 // not implemented yet
1505 // adding phase to the TreeS user info
1507 fLoader->TreeS()->GetUserInfo()->Add(new TParameter<float>("lhcphase0",lhcph));
1509 fTPCParam->SetZeroSup(0);
1511 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1512 if (IsSectorActive(isec)) {
1513 Hits2DigitsSector(isec);
1516 fLoader->WriteSDigits("OVERWRITE");
1518 //this line prevents the crash in the similar one
1519 //on the beginning of this method
1520 //destructor attempts to reset the tree, which is deleted by the loader
1521 //need to be redesign
1522 if(GetDigitsArray()) delete GetDigitsArray();
1523 SetDigitsArray(0x0);
1525 //__________________________________________________________________
1527 void AliTPC::Hits2SDigits()
1530 //-----------------------------------------------------------
1531 // summable digits - 16 bit "ADC", no noise, no saturation
1532 //-----------------------------------------------------------
1534 if (!fTPCParam->IsGeoRead()){
1536 // read transformation matrices for gGeoManager
1538 fTPCParam->ReadGeoMatrices();
1541 fLoader->LoadHits("read");
1542 fLoader->LoadSDigits("recreate");
1543 AliRunLoader* runLoader = fLoader->GetRunLoader();
1545 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1546 runLoader->GetEvent(iEvent);
1549 Hits2SDigits2(iEvent);
1552 fLoader->UnloadHits();
1553 fLoader->UnloadSDigits();
1554 if (fDebugStreamer) {
1555 delete fDebugStreamer;
1559 //_____________________________________________________________________________
1561 void AliTPC::Hits2DigitsSector(Int_t isec)
1563 //-------------------------------------------------------------------
1564 // TPC conversion from hits to digits.
1565 //-------------------------------------------------------------------
1567 //-----------------------------------------------------------------
1568 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1569 //-----------------------------------------------------------------
1571 //-------------------------------------------------------
1572 // Get the access to the track hits
1573 //-------------------------------------------------------
1575 // check if the parameters are set - important if one calls this method
1576 // directly, not from the Hits2Digits
1578 if(fDefaults == 0) SetDefaults();
1580 TTree *tH = fLoader->TreeH(); // pointer to the hits tree
1582 AliFatal("Can not find TreeH in folder");
1586 Stat_t ntracks = tH->GetEntries();
1588 Int_t nrows =fTPCParam->GetNRow(isec);
1590 TObjArray **row=new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1591 for(Int_t j=0;j<nrows+2;j++) row[j]=0;
1593 MakeSector(isec,nrows,tH,ntracks,row);
1595 //--------------------------------------------------------
1596 // Digitize this sector, row by row
1597 // row[i] is the pointer to the TObjArray of TVectors,
1598 // each one containing electrons accepted on this
1599 // row, assigned into tracks
1600 //--------------------------------------------------------
1604 if (fDigitsArray->GetTree()==0) {
1605 AliFatal("Tree not set in fDigitsArray");
1608 for (i=0;i<nrows;i++){
1610 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1612 DigitizeRow(i,isec,row);
1614 fDigitsArray->StoreRow(isec,i);
1616 Int_t ndig = dig->GetDigitSize();
1619 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1622 fDigitsArray->ClearRow(isec,i);
1625 } // end of the sector digitization
1627 for(i=0;i<nrows+2;i++){
1632 delete [] row; // delete the array of pointers to TObjArray-s
1635 } // end of Hits2DigitsSector
1638 //_____________________________________________________________________________
1639 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1641 //-----------------------------------------------------------
1642 // Single row digitization, coupling from the neighbouring
1643 // rows taken into account
1644 //-----------------------------------------------------------
1646 //-----------------------------------------------------------------
1647 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1648 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1649 //-----------------------------------------------------------------
1651 Float_t zerosup = fTPCParam->GetZeroSup();
1652 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor();
1653 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
1654 AliTPCCalROC * gainROC = gainTPC->GetCalROC(isec); // pad gains per given sector
1655 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(isec); // noise per given sector
1658 fCurrentIndex[1]= isec;
1661 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1662 Int_t nofTbins = fTPCParam->GetMaxTBin();
1663 Int_t indexRange[4];
1665 // Integrated signal for this row
1666 // and a single track signal
1669 TMatrixF *m1 = new TMatrixF(0,nofPads,0,nofTbins); // integrated
1670 TMatrixF *m2 = new TMatrixF(0,nofPads,0,nofTbins); // single
1672 TMatrixF &total = *m1;
1674 // Array of pointers to the label-signal list
1676 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1677 Float_t **pList = new Float_t* [nofDigits];
1681 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1687 for (Int_t row= row1;row<=row2;row++){
1688 Int_t nTracks= rows[row]->GetEntries();
1689 for (i1=0;i1<nTracks;i1++){
1690 fCurrentIndex[2]= row;
1691 fCurrentIndex[3]=irow+1;
1693 m2->Zero(); // clear single track signal matrix
1694 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1695 GetList(trackLabel,nofPads,m2,indexRange,pList);
1697 else GetSignal(rows[row],i1,0,m1,indexRange);
1703 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1705 Float_t fzerosup = zerosup+0.5;
1706 for(Int_t it=0;it<nofTbins;it++){
1707 for(Int_t ip=0;ip<nofPads;ip++){
1709 Float_t q=total(ip,it);
1710 if(fDigitsSwitch == 0){
1711 Float_t gain = gainROC->GetValue(irow,ip); // get gain for given - pad-row pad
1712 Float_t noisePad = noiseROC->GetValue(irow,ip);
1715 q+=GetNoise()*noisePad;
1716 if(q <=fzerosup) continue; // do not fill zeros
1718 if(q >= fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat() - 1; // saturation
1723 if(q <= 0.) continue; // do not fill zeros
1724 if(q>2000.) q=2000.;
1730 // "real" signal or electronic noise (list = -1)?
1733 for(Int_t j1=0;j1<3;j1++){
1734 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1739 <A NAME="AliDigits"></A>
1740 using of AliDigits object
1743 dig->SetDigitFast((Short_t)q,it,ip);
1744 if (fDigitsArray->IsSimulated()) {
1745 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1746 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1747 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1750 } // end of loop over time buckets
1751 } // end of lop over pads
1757 // glitch filters if normal simulated digits
1759 if(!fDigitsSwitch) ((AliSimDigits*)dig)->GlitchFilter();
1761 // This row has been digitized, delete nonused stuff
1764 for(lp=0;lp<nofDigits;lp++){
1765 if(pList[lp]) delete [] pList[lp];
1773 } // end of DigitizeRow
1775 //_____________________________________________________________________________
1777 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1778 TMatrixF *m1, TMatrixF *m2,Int_t *indexRange)
1781 //---------------------------------------------------------------
1782 // Calculates 2-D signal (pad,time) for a single track,
1783 // returns a pointer to the signal matrix and the track label
1784 // No digitization is performed at this level!!!
1785 //---------------------------------------------------------------
1787 //-----------------------------------------------------------------
1788 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1789 // Modified: Marian Ivanov
1790 //-----------------------------------------------------------------
1794 tv = (TVector*)p1->At(ntr); // pointer to a track
1797 Float_t label = v(0);
1798 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1))/2;
1800 Int_t nElectrons = (tv->GetNrows()-1)/5;
1801 indexRange[0]=9999; // min pad
1802 indexRange[1]=-1; // max pad
1803 indexRange[2]=9999; //min time
1804 indexRange[3]=-1; // max time
1806 TMatrixF &signal = *m1;
1807 TMatrixF &total = *m2;
1809 // Get LHC clock phase
1811 TParameter<float> *ph;
1812 if(fDigitsSwitch){// s-digits
1813 ph = (TParameter<float>*)fLoader->TreeS()->GetUserInfo()->FindObject("lhcphase0");
1815 else{ // normal digits
1816 ph = (TParameter<float>*)fLoader->TreeD()->GetUserInfo()->FindObject("lhcphase0");
1818 // Loop over all electrons
1820 for(Int_t nel=0; nel<nElectrons; nel++){
1822 Float_t aval = v(idx+4);
1823 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1824 Float_t xyz[4]={v(idx+1),v(idx+2),v(idx+3),v(idx+5)};
1825 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,
1826 fCurrentIndex[3],ph->GetVal());
1828 Int_t *index = fTPCParam->GetResBin(0);
1829 Float_t *weight = & (fTPCParam->GetResWeight(0));
1831 if (n>0) for (Int_t i =0; i<n; i++){
1832 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1835 Int_t time=index[2];
1836 Float_t qweight = *(weight)*eltoadcfac;
1838 if (m1!=0) signal(pad,time)+=qweight;
1839 total(pad,time)+=qweight;
1840 if (indexRange[0]>pad) indexRange[0]=pad;
1841 if (indexRange[1]<pad) indexRange[1]=pad;
1842 if (indexRange[2]>time) indexRange[2]=time;
1843 if (indexRange[3]<time) indexRange[3]=time;
1850 } // end of loop over electrons
1852 return label; // returns track label when finished
1855 //_____________________________________________________________________________
1856 void AliTPC::GetList(Float_t label,Int_t np,TMatrixF *m,
1857 Int_t *indexRange, Float_t **pList)
1859 //----------------------------------------------------------------------
1860 // Updates the list of tracks contributing to digits for a given row
1861 //----------------------------------------------------------------------
1863 //-----------------------------------------------------------------
1864 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1865 //-----------------------------------------------------------------
1867 TMatrixF &signal = *m;
1869 // lop over nonzero digits
1871 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1872 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1875 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1877 if(signal(ip,it)<0.5) continue;
1879 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1881 if(!pList[globalIndex]){
1884 // Create new list (6 elements - 3 signals and 3 labels),
1887 pList[globalIndex] = new Float_t [6];
1891 *pList[globalIndex] = -1.;
1892 *(pList[globalIndex]+1) = -1.;
1893 *(pList[globalIndex]+2) = -1.;
1894 *(pList[globalIndex]+3) = -1.;
1895 *(pList[globalIndex]+4) = -1.;
1896 *(pList[globalIndex]+5) = -1.;
1898 *pList[globalIndex] = label;
1899 *(pList[globalIndex]+3) = signal(ip,it);
1903 // check the signal magnitude
1905 Float_t highest = *(pList[globalIndex]+3);
1906 Float_t middle = *(pList[globalIndex]+4);
1907 Float_t lowest = *(pList[globalIndex]+5);
1910 // compare the new signal with already existing list
1913 if(signal(ip,it)<lowest) continue; // neglect this track
1917 if (signal(ip,it)>highest){
1918 *(pList[globalIndex]+5) = middle;
1919 *(pList[globalIndex]+4) = highest;
1920 *(pList[globalIndex]+3) = signal(ip,it);
1922 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1923 *(pList[globalIndex]+1) = *pList[globalIndex];
1924 *pList[globalIndex] = label;
1926 else if (signal(ip,it)>middle){
1927 *(pList[globalIndex]+5) = middle;
1928 *(pList[globalIndex]+4) = signal(ip,it);
1930 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1931 *(pList[globalIndex]+1) = label;
1934 *(pList[globalIndex]+5) = signal(ip,it);
1935 *(pList[globalIndex]+2) = label;
1939 } // end of loop over pads
1940 } // end of loop over time bins
1943 //___________________________________________________________________
1944 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1945 Stat_t ntracks,TObjArray **row)
1948 //-----------------------------------------------------------------
1949 // Prepares the sector digitization, creates the vectors of
1950 // tracks for each row of this sector. The track vector
1951 // contains the track label and the position of electrons.
1952 //-----------------------------------------------------------------
1955 // The trasport of the electrons through TPC drift volume
1956 // Drift (drift velocity + velocity map(not yet implemented)))
1957 // Application of the random processes (diffusion, gas gain)
1958 // Systematic effects (ExB effect in drift volume + ROCs)
1961 // Loop over primary electrons:
1962 // Creation of the secondary electrons
1963 // Loop over electrons (primary+ secondaries)
1964 // Global coordinate frame:
1965 // 1. Skip electrons if attached
1966 // 2. ExB effect in drift volume
1967 // a.) Simulation calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1968 // b.) Reconstruction - calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
1969 // 3. Generation of gas gain (Random - Exponential distribution)
1970 // 4. TransportElectron function (diffusion)
1972 // 5. Conversion to the local coordinate frame pad-row, pad, timebin
1973 // 6. Apply Time0 shift - AliTPCCalPad class
1974 // a.) Plus sign in simulation
1975 // b.) Minus sign in reconstruction
1978 //-----------------------------------------------------------------
1979 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1980 // Origin: Marian Ivanov, marian.ivanov@cern.ch
1981 //-----------------------------------------------------------------
1982 AliTPCcalibDB* const calib=AliTPCcalibDB::Instance();
1983 if (gAlice){ // Set correctly the magnetic field in the ExB calculation
1984 if (!calib->GetExB()){
1985 AliMagF * field = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField());
1987 calib->SetExBField(field);
1992 Float_t gasgain = fTPCParam->GetGasGain();
1993 gasgain = gasgain/fGainFactor;
1997 AliTPChit *tpcHit; // pointer to a sigle TPC hit
2000 if (fHitType>1) branch = TH->GetBranch("TPC2");
2001 else branch = TH->GetBranch("TPC");
2004 //----------------------------------------------
2005 // Create TObjArray-s, one for each row,
2006 // each TObjArray will store the TVectors
2007 // of electrons, one TVectors per each track.
2008 //----------------------------------------------
2010 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
2011 TVector **tracks = new TVector* [nrows+2]; //pointers to the track vectors
2013 for(i=0; i<nrows+2; i++){
2014 row[i] = new TObjArray;
2021 //--------------------------------------------------------------------
2022 // Loop over tracks, the "track" contains the full history
2023 //--------------------------------------------------------------------
2025 Int_t previousTrack,currentTrack;
2026 previousTrack = -1; // nothing to store so far!
2028 for(Int_t track=0;track<ntracks;track++){
2029 Bool_t isInSector=kTRUE;
2031 isInSector = TrackInVolume(isec,track);
2032 if (!isInSector) continue;
2034 branch->GetEntry(track); // get next track
2038 tpcHit = (AliTPChit*)FirstHit(-1);
2040 //--------------------------------------------------------------
2042 //--------------------------------------------------------------
2047 Int_t sector=tpcHit->fSector; // sector number
2049 tpcHit = (AliTPChit*) NextHit();
2053 // Remove hits which arrive before the TPC opening gate signal
2054 if(((fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
2055 /fTPCParam->GetDriftV()+tpcHit->Time())<fTPCParam->GetGateDelay()) {
2056 tpcHit = (AliTPChit*) NextHit();
2060 currentTrack = tpcHit->Track(); // track number
2062 if(currentTrack != previousTrack){
2064 // store already filled fTrack
2066 for(i=0;i<nrows+2;i++){
2067 if(previousTrack != -1){
2068 if(nofElectrons[i]>0){
2069 TVector &v = *tracks[i];
2070 v(0) = previousTrack;
2071 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
2072 row[i]->Add(tracks[i]);
2075 delete tracks[i]; // delete empty TVector
2081 tracks[i] = new TVector(601); // TVectors for the next fTrack
2083 } // end of loop over rows
2085 previousTrack=currentTrack; // update track label
2088 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2090 //---------------------------------------------------
2091 // Calculate the electron attachment probability
2092 //---------------------------------------------------
2095 Float_t time = 1.e6*(fTPCParam->GetZLength(isec)-TMath::Abs(tpcHit->Z()))
2096 /fTPCParam->GetDriftV();
2098 Float_t attProb = fTPCParam->GetAttCoef()*
2099 fTPCParam->GetOxyCont()*time; // fraction!
2101 //-----------------------------------------------
2102 // Loop over electrons
2103 //-----------------------------------------------
2106 for(Int_t nel=0;nel<qI;nel++){
2107 // skip if electron lost due to the attachment
2108 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2113 Double_t dxyz0[3],dxyz1[3];
2114 dxyz0[0]=tpcHit->X();
2115 dxyz0[1]=tpcHit->Y();
2116 dxyz0[2]=tpcHit->Z();
2117 if (calib->GetExB()){
2118 calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
2120 AliError("Not valid ExB calibration");
2121 dxyz1[0]=tpcHit->X();
2122 dxyz1[1]=tpcHit->Y();
2123 dxyz1[2]=tpcHit->Z();
2131 // protection for the nonphysical avalanche size (10**6 maximum)
2133 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2134 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2137 TransportElectron(xyz,index);
2139 Int_t padrow = fTPCParam->GetPadRow(xyz,index);
2141 // Add Time0 correction due unisochronity
2142 // xyz[0] - pad row coordinate
2143 // xyz[1] - pad coordinate
2144 // xyz[2] - is in now time bin coordinate system
2145 Float_t correction =0;
2146 if (calib->GetPadTime0()){
2147 if (!calib->GetPadTime0()->GetCalROC(isec)) continue;
2148 Int_t npads = fTPCParam->GetNPads(isec,padrow);
2149 // Int_t pad = TMath::Nint(xyz[1]+fTPCParam->GetNPads(isec,TMath::Nint(xyz[0]))*0.5);
2150 // pad numbering from -npads/2 .. npads/2-1
2151 Int_t pad = TMath::Nint(xyz[1]+npads/2);
2153 if (pad>=npads) pad=npads-1;
2154 correction = calib->GetPadTime0()->GetCalROC(isec)->GetValue(padrow,pad);
2155 // printf("%d\t%d\t%d\t%f\n",isec,padrow,pad,correction);
2156 if (fDebugStreamer){
2157 (*fDebugStreamer)<<"Time0"<<
2165 "cor="<<correction<<
2170 xyz[2]+=fTPCParam->GetNTBinsL1(); // adding Level 1 time bin offset
2172 // Electron track time (for pileup simulation)
2173 xyz[2]+=tpcHit->Time()/fTPCParam->GetTSample(); // adding time of flight
2177 // row 0 - cross talk from the innermost row
2178 // row fNRow+1 cross talk from the outermost row
2179 rowNumber = index[2]+1;
2180 //transform position to local digit coordinates
2181 //relative to nearest pad row
2182 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2184 if (isec <fTPCParam->GetNInnerSector()) {
2185 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2186 y1 = fTPCParam->GetYInner(rowNumber);
2189 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2190 y1 = fTPCParam->GetYOuter(rowNumber);
2192 // gain inefficiency at the wires edges - linear
2195 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.)); */
2197 nofElectrons[rowNumber]++;
2198 //----------------------------------
2199 // Expand vector if necessary
2200 //----------------------------------
2201 if(nofElectrons[rowNumber]>120){
2202 Int_t range = tracks[rowNumber]->GetNrows();
2203 if((nofElectrons[rowNumber])>(range-1)/5){
2205 tracks[rowNumber]->ResizeTo(range+500); // Add 100 electrons
2209 TVector &v = *tracks[rowNumber];
2210 Int_t idx = 5*nofElectrons[rowNumber]-4;
2211 Real_t * position = &(((TVector&)v)(idx)); //make code faster
2212 memcpy(position,xyz,5*sizeof(Float_t));
2214 } // end of loop over electrons
2216 tpcHit = (AliTPChit*)NextHit();
2218 } // end of loop over hits
2219 } // end of loop over tracks
2222 // store remaining track (the last one) if not empty
2225 for(i=0;i<nrows+2;i++){
2226 if(nofElectrons[i]>0){
2227 TVector &v = *tracks[i];
2228 v(0) = previousTrack;
2229 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
2230 row[i]->Add(tracks[i]);
2239 delete [] nofElectrons;
2241 } // end of MakeSector
2244 //_____________________________________________________________________________
2248 // Initialise TPC detector after definition of geometry
2250 AliDebug(1,"*********************************************");
2253 //_____________________________________________________________________________
2254 void AliTPC::ResetDigits()
2257 // Reset number of digits and the digits array for this detector
2260 if (fDigits) fDigits->Clear();
2265 //_____________________________________________________________________________
2266 void AliTPC::SetSens(Int_t sens)
2269 //-------------------------------------------------------------
2270 // Activates/deactivates the sensitive strips at the center of
2271 // the pad row -- this is for the space-point resolution calculations
2272 //-------------------------------------------------------------
2274 //-----------------------------------------------------------------
2275 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2276 //-----------------------------------------------------------------
2282 void AliTPC::SetSide(Float_t side=0.)
2284 // choice of the TPC side
2289 //_____________________________________________________________________________
2291 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2294 // electron transport taking into account:
2296 // 2.ExB at the wires
2297 // 3. nonisochronity
2299 // xyz and index must be already transformed to system 1
2302 fTPCParam->Transform1to2(xyz,index); // mis-alignment applied in this step
2305 Float_t driftl=xyz[2];
2306 if(driftl<0.01) driftl=0.01;
2307 driftl=TMath::Sqrt(driftl);
2308 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2309 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2310 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2311 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2312 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2316 if (fTPCParam->GetMWPCReadout()==kTRUE){
2317 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2318 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2320 //add nonisochronity (not implemented yet)
2326 //______________________________________________________________________
2327 AliTPChit::AliTPChit()
2339 //_____________________________________________________________________________
2340 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits)
2341 :AliHit(shunt,track),
2348 // Creates a TPC hit object
2359 //________________________________________________________________________
2360 // Additional code because of the AliTPCTrackHitsV2
2362 void AliTPC::MakeBranch(Option_t *option)
2365 // Create a new branch in the current Root Tree
2366 // The branch of fHits is automatically split
2367 // MI change 14.09.2000
2369 if (fHitType<2) return;
2370 char branchname[10];
2371 //sprintf(branchname,"%s2",GetName());
2372 snprintf(branchname,10,"%s2",GetName());
2374 // Get the pointer to the header
2375 const char *cH = strstr(option,"H");
2377 if (fTrackHits && fLoader->TreeH() && cH && fHitType&4) {
2378 AliDebug(1,"Making branch for Type 4 Hits");
2379 fLoader->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2382 // if (fTrackHitsOld && fLoader->TreeH() && cH && fHitType&2) {
2383 // AliDebug(1,"Making branch for Type 2 Hits");
2384 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2385 // fLoader->TreeH(),fBufferSize,99);
2386 // fLoader->TreeH()->GetListOfBranches()->Add(branch);
2390 void AliTPC::SetTreeAddress()
2392 //Sets tree address for hits
2394 if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2395 AliDetector::SetTreeAddress();
2397 if (fHitType>1) SetTreeAddress2();
2400 void AliTPC::SetTreeAddress2()
2403 // Set branch address for the TrackHits Tree
2408 char branchname[20];
2409 //sprintf(branchname,"%s2",GetName());
2410 snprintf(branchname,20,"%s2",GetName());
2412 // Branch address for hit tree
2413 TTree *treeH = fLoader->TreeH();
2414 if ((treeH)&&(fHitType&4)) {
2415 branch = treeH->GetBranch(branchname);
2417 branch->SetAddress(&fTrackHits);
2418 AliDebug(1,"fHitType&4 Setting");
2421 AliDebug(1,"fHitType&4 Failed (can not find branch)");
2424 // if ((treeH)&&(fHitType&2)) {
2425 // branch = treeH->GetBranch(branchname);
2427 // branch->SetAddress(&fTrackHitsOld);
2428 // AliDebug(1,"fHitType&2 Setting");
2431 // AliDebug(1,"fHitType&2 Failed (can not find branch)");
2435 void AliTPC::FinishPrimary()
2437 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2438 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2442 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2445 // add hit to the list
2449 int primary = gAlice->GetMCApp()->GetPrimary(track);
2450 gAlice->GetMCApp()->Particle(primary)->SetBit(kKeepBit);
2454 gAlice->GetMCApp()->FlagTrack(track);
2456 if (fTrackHits && fHitType&4)
2457 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2458 hits[1],hits[2],(Int_t)hits[3],hits[4]);
2459 // if (fTrackHitsOld &&fHitType&2 )
2460 // fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2461 // hits[1],hits[2],(Int_t)hits[3]);
2465 void AliTPC::ResetHits()
2467 if (fHitType&1) AliDetector::ResetHits();
2468 if (fHitType>1) ResetHits2();
2471 void AliTPC::ResetHits2()
2475 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2476 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2480 AliHit* AliTPC::FirstHit(Int_t track)
2482 if (fHitType>1) return FirstHit2(track);
2483 return AliDetector::FirstHit(track);
2485 AliHit* AliTPC::NextHit()
2490 if (fHitType>1) return NextHit2();
2492 return AliDetector::NextHit();
2495 AliHit* AliTPC::FirstHit2(Int_t track)
2498 // Initialise the hit iterator
2499 // Return the address of the first hit for track
2500 // If track>=0 the track is read from disk
2501 // while if track<0 the first hit of the current
2502 // track is returned
2505 gAlice->GetMCApp()->ResetHits();
2506 fLoader->TreeH()->GetEvent(track);
2509 if (fTrackHits && fHitType&4) {
2510 fTrackHits->First();
2511 return fTrackHits->GetHit();
2513 // if (fTrackHitsOld && fHitType&2) {
2514 // fTrackHitsOld->First();
2515 // return fTrackHitsOld->GetHit();
2521 AliHit* AliTPC::NextHit2()
2524 //Return the next hit for the current track
2527 // if (fTrackHitsOld && fHitType&2) {
2528 // fTrackHitsOld->Next();
2529 // return fTrackHitsOld->GetHit();
2533 return fTrackHits->GetHit();
2539 void AliTPC::RemapTrackHitIDs(Int_t *map)
2544 if (!fTrackHits) return;
2546 // if (fTrackHitsOld && fHitType&2){
2547 // AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2548 // for (UInt_t i=0;i<arr->GetSize();i++){
2549 // AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2550 // info->fTrackID = map[info->fTrackID];
2553 // if (fTrackHitsOld && fHitType&4){
2554 if (fTrackHits && fHitType&4){
2555 TClonesArray * arr = fTrackHits->GetArray();;
2556 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2557 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2558 info->SetTrackID(map[info->GetTrackID()]);
2563 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2565 //return bool information - is track in given volume
2566 //load only part of the track information
2567 //return true if current track is in volume
2570 // if (fTrackHitsOld && fHitType&2) {
2571 // TBranch * br = fLoader->TreeH()->GetBranch("fTrackHitsInfo");
2572 // br->GetEvent(track);
2573 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2574 // for (UInt_t j=0;j<ar->GetSize();j++){
2575 // if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2579 if (fTrackHits && fHitType&4) {
2580 TBranch * br1 = fLoader->TreeH()->GetBranch("fVolumes");
2581 TBranch * br2 = fLoader->TreeH()->GetBranch("fNVolumes");
2582 br2->GetEvent(track);
2583 br1->GetEvent(track);
2584 Int_t *volumes = fTrackHits->GetVolumes();
2585 Int_t nvolumes = fTrackHits->GetNVolumes();
2586 if (!volumes && nvolumes>0) {
2587 AliWarning(Form("Problematic track\t%d\t%d",track,nvolumes));
2590 for (Int_t j=0;j<nvolumes; j++)
2591 if (volumes[j]==id) return kTRUE;
2595 TBranch * br = fLoader->TreeH()->GetBranch("fSector");
2596 br->GetEvent(track);
2597 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2598 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2606 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
2609 fLoader = new AliTPCLoader(GetName(),topfoldername);
2613 ////////////////////////////////////////////////////////////////////////
2614 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
2616 // load TPC paarmeters from a given file or create new if the object
2617 // is not found there
2618 // 12/05/2003 This method should be moved to the AliTPCLoader
2619 // and one has to decide where to store the TPC parameters
2622 //sprintf(paramName,"75x40_100x60_150x60");
2623 snprintf(paramName,50,"75x40_100x60_150x60");
2624 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
2626 AliDebugClass(1,Form("TPC parameters %s found.",paramName));
2628 AliWarningClass("TPC parameters not found. Create new (they may be incorrect)");
2629 //paramTPC = new AliTPCParamSR;
2630 paramTPC = AliTPCcalibDB::Instance()->GetParameters();
2631 if (!paramTPC->IsGeoRead()){
2633 // read transformation matrices for gGeoManager
2635 paramTPC->ReadGeoMatrices();
2641 // the older version of parameters can be accessed with this code.
2642 // In some cases, we have old parameters saved in the file but
2643 // digits were created with new parameters, it can be distinguish
2644 // by the name of TPC TreeD. The code here is just for the case
2645 // we would need to compare with old data, uncomment it if needed.
2647 // char paramName[50];
2648 // sprintf(paramName,"75x40_100x60");
2649 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
2651 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2653 // sprintf(paramName,"75x40_100x60_150x60");
2654 // paramTPC=(AliTPCParam*)in->Get(paramName);
2656 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2658 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
2660 // paramTPC = new AliTPCParamSR;