1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // Time Projection Chamber //
21 // This class contains the basic functions for the Time Projection Chamber //
22 // detector. Functions specific to one particular geometry are //
23 // contained in the derived classes //
27 <img src="picts/AliTPCClass.gif">
32 ///////////////////////////////////////////////////////////////////////////////
36 #include <Riostream.h>
40 #include <TGeometry.h>
41 #include <TInterpreter.h>
46 #include <TObjectTable.h>
47 #include <TParticle.h>
53 #include <TVirtualMC.h>
56 #include <TStopwatch.h>
58 #include "AliDigits.h"
60 #include "AliPoints.h"
62 #include "AliRunLoader.h"
63 #include "AliSimDigits.h"
66 #include "AliTPCDigitsArray.h"
67 #include "AliTPCLoader.h"
68 #include "AliTPCPRF2D.h"
69 #include "AliTPCParamSR.h"
70 #include "AliTPCRF1D.h"
71 //#include "AliTPCTrackHits.h"
72 #include "AliTPCTrackHitsV2.h"
73 #include "AliTrackReference.h"
75 #include "AliTPCDigitizer.h"
76 #include "AliTPCBuffer.h"
77 #include "AliTPCDDLRawData.h"
82 //_____________________________________________________________________________
86 // Default constructor
96 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
97 fHitType = 4; // ROOT containers
99 fHitType = 2; //default CONTAINERS - based on ROOT structure
107 //_____________________________________________________________________________
108 AliTPC::AliTPC(const char *name, const char *title)
109 : AliDetector(name,title)
112 // Standard constructor
116 // Initialise arrays of hits and digits
117 fHits = new TClonesArray("AliTPChit", 176);
118 gAlice->GetMCApp()->AddHitList(fHits);
122 fTrackHits = new AliTPCTrackHitsV2;
123 fTrackHits->SetHitPrecision(0.002);
124 fTrackHits->SetStepPrecision(0.003);
125 fTrackHits->SetMaxDistance(100);
127 //fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
128 //fTrackHitsOld->SetHitPrecision(0.002);
129 //fTrackHitsOld->SetStepPrecision(0.003);
130 //fTrackHitsOld->SetMaxDistance(100);
134 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
135 fHitType = 4; // ROOT containers
141 // Initialise counters
147 // Initialise color attributes
148 SetMarkerColor(kYellow);
151 // Set TPC parameters
155 if (!strcmp(title,"Default")) {
156 fTPCParam = new AliTPCParamSR;
158 AliWarning("In Config.C you must set non-default parameters.");
164 //_____________________________________________________________________________
165 AliTPC::AliTPC(const AliTPC& t):AliDetector(t){
167 // dummy copy constructor
180 delete fTrackHits; //MI 15.09.2000
181 // delete fTrackHitsOld; //MI 10.12.2001
182 if (fNoiseTable) delete [] fNoiseTable;
186 //_____________________________________________________________________________
187 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
190 // Add a hit to the list
193 TClonesArray &lhits = *fHits;
194 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
197 AddHit2(track,vol,hits);
200 //_____________________________________________________________________________
201 void AliTPC::BuildGeometry()
205 // Build TPC ROOT TNode geometry for the event display
210 const int kColorTPC=19;
211 char name[5], title[25];
212 const Double_t kDegrad=TMath::Pi()/180;
213 const Double_t kRaddeg=180./TMath::Pi();
216 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
217 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
219 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
220 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
222 Int_t nLo = fTPCParam->GetNInnerSector()/2;
223 Int_t nHi = fTPCParam->GetNOuterSector()/2;
225 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
226 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
227 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
228 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
231 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
232 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
238 // Get ALICE top node
241 nTop=gAlice->GetGeometry()->GetNode("alice");
245 rl = fTPCParam->GetInnerRadiusLow();
246 ru = fTPCParam->GetInnerRadiusUp();
250 sprintf(name,"LS%2.2d",i);
252 sprintf(title,"TPC low sector %3d",i);
255 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
256 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
257 tubs->SetNumberOfDivisions(1);
259 nNode = new TNode(name,title,name,0,0,0,"");
260 nNode->SetLineColor(kColorTPC);
266 rl = fTPCParam->GetOuterRadiusLow();
267 ru = fTPCParam->GetOuterRadiusUp();
270 sprintf(name,"US%2.2d",i);
272 sprintf(title,"TPC upper sector %d",i);
274 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
275 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
276 tubs->SetNumberOfDivisions(1);
278 nNode = new TNode(name,title,name,0,0,0,"");
279 nNode->SetLineColor(kColorTPC);
285 //_____________________________________________________________________________
286 void AliTPC::CreateMaterials()
288 //-----------------------------------------------
289 // Create Materials for for TPC simulations
290 //-----------------------------------------------
292 //-----------------------------------------------------------------
293 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
294 //-----------------------------------------------------------------
296 Int_t iSXFLD=gAlice->Field()->Integ();
297 Float_t sXMGMX=gAlice->Field()->Max();
299 Float_t amat[5]; // atomic numbers
300 Float_t zmat[5]; // z
301 Float_t wmat[5]; // proportions
307 //***************** Gases *************************
309 //-------------------------------------------------
311 //-------------------------------------------------
322 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
332 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
335 //--------------------------------------------------------------
337 //--------------------------------------------------------------
354 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
356 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
371 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
373 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
389 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
391 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
393 //----------------------------------------------------------------
394 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
395 //----------------------------------------------------------------
401 Float_t rho,absl,x0,buf[1];
405 for(nc = 0;nc<fNoComp;nc++) {
407 // retrive material constants
409 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,x0,absl,buf,nbuf);
414 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
416 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
417 density += fMixtProp[nc]*rho; // density of the mixture
421 // mixture proportions by weight!
423 for(nc = 0;nc<fNoComp;nc++) {
425 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
427 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
428 apure[nnc] : amol[nnc])/am;
432 // Drift gases 1 - nonsensitive, 2 - sensitive
434 // AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
435 // AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
445 wmat[2] = wmat[1]*0.728;
450 AliMixture(31,"Drift gas 1",amat,zmat,density,3,wmat);
451 AliMixture(32,"Drift gas 2",amat,zmat,density,3,wmat);
460 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
476 AliMixture(24,"Air",amat,zmat,density,2,wmat);
478 //----------------------------------------------------------------------
480 //----------------------------------------------------------------------
502 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
524 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
542 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
560 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
562 // SiO2 - used later for the glass fiber
574 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
583 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
592 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
601 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
619 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
638 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
640 // Epoxy - C14 H20 O3
657 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
665 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
669 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,x0,absl,buf,nbuf);
683 wmat[0]=0.644; // by weight! C
684 wmat[1]=0.356; // epoxy altogether
687 wmat[3]=wmat[1]*0.203;
688 wmat[2]=wmat[1]*0.085;
692 density=0.5*(1.25+2.265);
694 AliMixture(47,"Cfiber",amat,zmat,density,4,wmat);
698 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,x0,absl,buf,nbuf);
718 wmat[0]=0.725; // by weight!
719 wmat[1]=wmat[0]*0.533;
723 wmat[3]=wmat[2]*0.085;
724 wmat[4]=wmat[2]*0.203;
729 AliMixture(39,"G10",amat,zmat,density,5,wmat);
734 //----------------------------------------------------------
735 // tracking media for gases
736 //----------------------------------------------------------
738 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
739 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
740 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
741 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
743 //-----------------------------------------------------------
744 // tracking media for solids
745 //-----------------------------------------------------------
747 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
748 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
749 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
750 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
751 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
752 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
753 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
754 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
755 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
756 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
757 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
758 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
762 void AliTPC::GenerNoise(Int_t tablesize)
765 //Generate table with noise
772 if (fNoiseTable) delete[] fNoiseTable;
773 fNoiseTable = new Float_t[tablesize];
774 fNoiseDepth = tablesize;
775 fCurrentNoise =0; //!index of the noise in the noise table
777 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
778 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
781 Float_t AliTPC::GetNoise()
783 // get noise from table
784 // if ((fCurrentNoise%10)==0)
785 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
786 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
787 return fNoiseTable[fCurrentNoise++];
788 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
792 Bool_t AliTPC::IsSectorActive(Int_t sec) const
795 // check if the sector is active
796 if (!fActiveSectors) return kTRUE;
797 else return fActiveSectors[sec];
800 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
802 // activate interesting sectors
803 SetTreeAddress();//just for security
804 if (fActiveSectors) delete [] fActiveSectors;
805 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
806 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
807 for (Int_t i=0;i<n;i++)
808 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
812 void AliTPC::SetActiveSectors(Int_t flag)
815 // activate sectors which were hitted by tracks
817 SetTreeAddress();//just for security
818 if (fHitType==0) return; // if Clones hit - not short volume ID information
819 if (fActiveSectors) delete [] fActiveSectors;
820 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
822 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
825 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
829 AliFatal("Can not find TreeH in folder");
832 if (fHitType>1) branch = TreeH()->GetBranch("TPC2");
833 else branch = TreeH()->GetBranch("TPC");
834 Stat_t ntracks = TreeH()->GetEntries();
835 // loop over all hits
836 AliDebug(1,Form("Got %d tracks",ntracks));
838 for(Int_t track=0;track<ntracks;track++) {
841 if (fTrackHits && fHitType&4) {
842 TBranch * br1 = TreeH()->GetBranch("fVolumes");
843 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
844 br1->GetEvent(track);
845 br2->GetEvent(track);
846 Int_t *volumes = fTrackHits->GetVolumes();
847 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
848 fActiveSectors[volumes[j]]=kTRUE;
852 // if (fTrackHitsOld && fHitType&2) {
853 // TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
854 // br->GetEvent(track);
855 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
856 // for (UInt_t j=0;j<ar->GetSize();j++){
857 // fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
866 //_____________________________________________________________________________
867 void AliTPC::Digits2Raw()
869 // convert digits of the current event to raw data
871 static const Int_t kThreshold = 0;
872 static const Bool_t kCompress = kTRUE;
874 fLoader->LoadDigits();
875 TTree* digits = fLoader->TreeD();
877 AliError("No digits tree");
882 AliSimDigits* digrow = &digarr;
883 digits->GetBranch("Segment")->SetAddress(&digrow);
885 const char* fileName = "AliTPCDDL.dat";
886 AliTPCBuffer* buffer = new AliTPCBuffer(fileName);
890 // 2: txt files with digits
891 //BE CAREFUL, verbose level 2 MUST be used only for debugging and
892 //it is highly suggested to use this mode only for debugging digits files
893 //reasonably small, because otherwise the size of the txt files can reach
894 //quickly several MB wasting time and disk space.
895 buffer->SetVerbose(0);
897 Int_t nEntries = Int_t(digits->GetEntries());
898 Int_t previousSector = -1;
900 for (Int_t i = 0; i < nEntries; i++) {
903 fTPCParam->AdjustSectorRow(digarr.GetID(), sector, row);
904 if(previousSector != sector) {
906 previousSector = sector;
909 if (sector < 36) { //inner sector [0;35]
911 //the whole row is written into the output file
912 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
913 sector, subSector, row);
915 //only the pads in the range [37;48] are written into the output file
916 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 1,
917 sector, subSector, row);
919 //only the pads outside the range [37;48] are written into the output file
920 buffer->WriteRowBinary(kThreshold, digrow, 37, 48, 2,
921 sector, subSector, row);
924 } else { //outer sector [36;71]
925 if (row == 54) subSector = 2;
926 if ((row != 27) && (row != 76)) {
927 buffer->WriteRowBinary(kThreshold, digrow, 0, 0, 0,
928 sector, subSector, row);
929 } else if (row == 27) {
930 //only the pads outside the range [43;46] are written into the output file
931 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 2,
932 sector, subSector, row);
934 //only the pads in the range [43;46] are written into the output file
935 buffer->WriteRowBinary(kThreshold, digrow, 43, 46, 1,
936 sector, subSector, row);
937 } else if (row == 76) {
938 //only the pads outside the range [33;88] are written into the output file
939 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 2,
940 sector, subSector, row);
942 //only the pads in the range [33;88] are written into the output file
943 buffer->WriteRowBinary(kThreshold, digrow, 33, 88, 1,
944 sector, subSector, row);
950 fLoader->UnloadDigits();
952 AliTPCDDLRawData rawWriter;
953 rawWriter.SetVerbose(0);
955 rawWriter.RawData(fileName);
956 gSystem->Unlink(fileName);
959 AliInfo("Compressing raw data");
960 rawWriter.RawDataCompDecompress(kTRUE);
961 gSystem->Unlink("Statistics");
967 //______________________________________________________________________
968 AliDigitizer* AliTPC::CreateDigitizer(AliRunDigitizer* manager) const
970 return new AliTPCDigitizer(manager);
973 void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
975 //create digits from summable digits
976 GenerNoise(500000); //create teble with noise
978 //conect tree with sSDigits
979 TTree *t = fLoader->TreeS();
982 fLoader->LoadSDigits("READ");
983 t = fLoader->TreeS();
985 AliError("Can not get input TreeS");
990 if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
992 AliSimDigits digarr, *dummy=&digarr;
993 TBranch* sdb = t->GetBranch("Segment");
995 AliError("Can not find branch with segments in TreeS.");
999 sdb->SetAddress(&dummy);
1001 Stat_t nentries = t->GetEntries();
1003 // set zero suppression
1005 fTPCParam->SetZeroSup(2);
1007 // get zero suppression
1009 Int_t zerosup = fTPCParam->GetZeroSup();
1011 //make tree with digits
1013 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1014 arr->SetClass("AliSimDigits");
1015 arr->Setup(fTPCParam);
1016 arr->MakeTree(fLoader->TreeD());
1018 AliTPCParam * par = fTPCParam;
1020 //Loop over segments of the TPC
1022 for (Int_t n=0; n<nentries; n++) {
1025 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1026 AliWarning(Form("Invalid segment ID ! %d",digarr.GetID()));
1029 if (!IsSectorActive(sec)) continue;
1031 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1032 Int_t nrows = digrow->GetNRows();
1033 Int_t ncols = digrow->GetNCols();
1035 digrow->ExpandBuffer();
1036 digarr.ExpandBuffer();
1037 digrow->ExpandTrackBuffer();
1038 digarr.ExpandTrackBuffer();
1041 Short_t * pamp0 = digarr.GetDigits();
1042 Int_t * ptracks0 = digarr.GetTracks();
1043 Short_t * pamp1 = digrow->GetDigits();
1044 Int_t * ptracks1 = digrow->GetTracks();
1045 Int_t nelems =nrows*ncols;
1046 Int_t saturation = fTPCParam->GetADCSat();
1047 //use internal structure of the AliDigits - for speed reason
1048 //if you cahnge implementation
1049 //of the Alidigits - it must be rewriten -
1050 for (Int_t i= 0; i<nelems; i++){
1051 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1053 if (q>saturation) q=saturation;
1056 ptracks1[0]=ptracks0[0];
1057 ptracks1[nelems]=ptracks0[nelems];
1058 ptracks1[2*nelems]=ptracks0[2*nelems];
1066 arr->StoreRow(sec,row);
1067 arr->ClearRow(sec,row);
1072 fLoader->WriteDigits("OVERWRITE");
1076 //__________________________________________________________________
1077 void AliTPC::SetDefaults(){
1079 // setting the defaults
1082 // Set response functions
1085 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1087 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1089 AliInfo("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...");
1091 param = new AliTPCParamSR();
1094 param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1097 AliFatal("No TPC parameters found");
1101 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1102 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1103 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1104 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1105 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1106 rf->SetOffset(3*param->GetZSigma());
1109 TDirectory *savedir=gDirectory;
1110 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1112 AliFatal("Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !");
1115 prfinner->Read("prf_07504_Gati_056068_d02");
1116 //PH Set different names
1117 s=prfinner->GetGRF()->GetName();
1119 prfinner->GetGRF()->SetName(s.Data());
1121 prfouter1->Read("prf_10006_Gati_047051_d03");
1122 s=prfouter1->GetGRF()->GetName();
1124 prfouter1->GetGRF()->SetName(s.Data());
1126 prfouter2->Read("prf_15006_Gati_047051_d03");
1127 s=prfouter2->GetGRF()->GetName();
1129 prfouter2->GetGRF()->SetName(s.Data());
1134 param->SetInnerPRF(prfinner);
1135 param->SetOuter1PRF(prfouter1);
1136 param->SetOuter2PRF(prfouter2);
1137 param->SetTimeRF(rf);
1147 //__________________________________________________________________
1148 void AliTPC::Hits2Digits()
1151 // creates digits from hits
1154 fLoader->LoadHits("read");
1155 fLoader->LoadDigits("recreate");
1156 AliRunLoader* runLoader = fLoader->GetRunLoader();
1158 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1159 runLoader->GetEvent(iEvent);
1161 Hits2Digits(iEvent);
1164 fLoader->UnloadHits();
1165 fLoader->UnloadDigits();
1167 //__________________________________________________________________
1168 void AliTPC::Hits2Digits(Int_t eventnumber)
1170 //----------------------------------------------------
1171 // Loop over all sectors for a single event
1172 //----------------------------------------------------
1173 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::GetRunLoaderName());
1174 rl->GetEvent(eventnumber);
1175 if (fLoader->TreeH() == 0x0) {
1176 if(fLoader->LoadHits()) {
1177 AliError("Can not load hits.");
1182 if (fLoader->TreeD() == 0x0 ) {
1183 fLoader->MakeTree("D");
1184 if (fLoader->TreeD() == 0x0 ) {
1185 AliError("Can not get TreeD");
1190 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1191 GenerNoise(500000); //create teble with noise
1193 //setup TPCDigitsArray
1195 if(GetDigitsArray()) delete GetDigitsArray();
1197 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1198 arr->SetClass("AliSimDigits");
1199 arr->Setup(fTPCParam);
1201 arr->MakeTree(fLoader->TreeD());
1202 SetDigitsArray(arr);
1204 fDigitsSwitch=0; // standard digits
1206 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1207 if (IsSectorActive(isec)) {
1208 AliDebug(1,Form("Hits2Digits","Sector %d is active.",isec));
1209 Hits2DigitsSector(isec);
1212 AliDebug(1,Form("Hits2Digits","Sector %d is NOT active.",isec));
1215 fLoader->WriteDigits("OVERWRITE");
1217 //this line prevents the crash in the similar one
1218 //on the beginning of this method
1219 //destructor attempts to reset the tree, which is deleted by the loader
1220 //need to be redesign
1221 if(GetDigitsArray()) delete GetDigitsArray();
1222 SetDigitsArray(0x0);
1226 //__________________________________________________________________
1227 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1230 //-----------------------------------------------------------
1231 // summable digits - 16 bit "ADC", no noise, no saturation
1232 //-----------------------------------------------------------
1234 //----------------------------------------------------
1235 // Loop over all sectors for a single event
1236 //----------------------------------------------------
1238 AliRunLoader* rl = fLoader->GetRunLoader();
1240 rl->GetEvent(eventnumber);
1241 if (fLoader->TreeH() == 0x0) {
1242 if(fLoader->LoadHits()) {
1243 AliError("Can not load hits.");
1250 if (fLoader->TreeS() == 0x0 ) {
1251 fLoader->MakeTree("S");
1254 if(fDefaults == 0) SetDefaults();
1256 GenerNoise(500000); //create table with noise
1257 //setup TPCDigitsArray
1259 if(GetDigitsArray()) delete GetDigitsArray();
1262 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1263 arr->SetClass("AliSimDigits");
1264 arr->Setup(fTPCParam);
1265 arr->MakeTree(fLoader->TreeS());
1267 SetDigitsArray(arr);
1269 fDigitsSwitch=1; // summable digits
1271 // set zero suppression to "0"
1273 fTPCParam->SetZeroSup(0);
1275 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1276 if (IsSectorActive(isec)) {
1277 Hits2DigitsSector(isec);
1280 fLoader->WriteSDigits("OVERWRITE");
1282 //this line prevents the crash in the similar one
1283 //on the beginning of this method
1284 //destructor attempts to reset the tree, which is deleted by the loader
1285 //need to be redesign
1286 if(GetDigitsArray()) delete GetDigitsArray();
1287 SetDigitsArray(0x0);
1289 //__________________________________________________________________
1291 void AliTPC::Hits2SDigits()
1294 //-----------------------------------------------------------
1295 // summable digits - 16 bit "ADC", no noise, no saturation
1296 //-----------------------------------------------------------
1298 fLoader->LoadHits("read");
1299 fLoader->LoadSDigits("recreate");
1300 AliRunLoader* runLoader = fLoader->GetRunLoader();
1302 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1303 runLoader->GetEvent(iEvent);
1306 Hits2SDigits2(iEvent);
1309 fLoader->UnloadHits();
1310 fLoader->UnloadSDigits();
1312 //_____________________________________________________________________________
1314 void AliTPC::Hits2DigitsSector(Int_t isec)
1316 //-------------------------------------------------------------------
1317 // TPC conversion from hits to digits.
1318 //-------------------------------------------------------------------
1320 //-----------------------------------------------------------------
1321 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1322 //-----------------------------------------------------------------
1324 //-------------------------------------------------------
1325 // Get the access to the track hits
1326 //-------------------------------------------------------
1328 // check if the parameters are set - important if one calls this method
1329 // directly, not from the Hits2Digits
1331 if(fDefaults == 0) SetDefaults();
1333 TTree *tH = TreeH(); // pointer to the hits tree
1335 AliFatal("Can not find TreeH in folder");
1339 Stat_t ntracks = tH->GetEntries();
1343 //-------------------------------------------
1344 // Only if there are any tracks...
1345 //-------------------------------------------
1349 Int_t nrows =fTPCParam->GetNRow(isec);
1351 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1353 MakeSector(isec,nrows,tH,ntracks,row);
1355 //--------------------------------------------------------
1356 // Digitize this sector, row by row
1357 // row[i] is the pointer to the TObjArray of TVectors,
1358 // each one containing electrons accepted on this
1359 // row, assigned into tracks
1360 //--------------------------------------------------------
1364 if (fDigitsArray->GetTree()==0) {
1365 AliFatal("Tree not set in fDigitsArray");
1368 for (i=0;i<nrows;i++){
1370 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1372 DigitizeRow(i,isec,row);
1374 fDigitsArray->StoreRow(isec,i);
1376 Int_t ndig = dig->GetDigitSize();
1379 Form("*** Sector, row, compressed digits %d %d %d ***\n",
1382 fDigitsArray->ClearRow(isec,i);
1385 } // end of the sector digitization
1387 for(i=0;i<nrows+2;i++){
1392 delete [] row; // delete the array of pointers to TObjArray-s
1396 } // end of Hits2DigitsSector
1399 //_____________________________________________________________________________
1400 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1402 //-----------------------------------------------------------
1403 // Single row digitization, coupling from the neighbouring
1404 // rows taken into account
1405 //-----------------------------------------------------------
1407 //-----------------------------------------------------------------
1408 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1409 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1410 //-----------------------------------------------------------------
1412 Float_t zerosup = fTPCParam->GetZeroSup();
1414 fCurrentIndex[1]= isec;
1417 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1418 Int_t nofTbins = fTPCParam->GetMaxTBin();
1419 Int_t indexRange[4];
1421 // Integrated signal for this row
1422 // and a single track signal
1425 TMatrixF *m1 = new TMatrixF(0,nofPads,0,nofTbins); // integrated
1426 TMatrixF *m2 = new TMatrixF(0,nofPads,0,nofTbins); // single
1428 TMatrixF &total = *m1;
1430 // Array of pointers to the label-signal list
1432 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1433 Float_t **pList = new Float_t* [nofDigits];
1437 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1443 for (Int_t row= row1;row<=row2;row++){
1444 Int_t nTracks= rows[row]->GetEntries();
1445 for (i1=0;i1<nTracks;i1++){
1446 fCurrentIndex[2]= row;
1447 fCurrentIndex[3]=irow+1;
1449 m2->Zero(); // clear single track signal matrix
1450 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1451 GetList(trackLabel,nofPads,m2,indexRange,pList);
1453 else GetSignal(rows[row],i1,0,m1,indexRange);
1459 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1461 Float_t fzerosup = zerosup+0.5;
1462 for(Int_t it=0;it<nofTbins;it++){
1463 for(Int_t ip=0;ip<nofPads;ip++){
1465 Float_t q=total(ip,it);
1466 if(fDigitsSwitch == 0){
1468 if(q <=fzerosup) continue; // do not fill zeros
1470 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1475 if(q <= 0.) continue; // do not fill zeros
1476 if(q>2000.) q=2000.;
1482 // "real" signal or electronic noise (list = -1)?
1485 for(Int_t j1=0;j1<3;j1++){
1486 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1491 <A NAME="AliDigits"></A>
1492 using of AliDigits object
1495 dig->SetDigitFast((Short_t)q,it,ip);
1496 if (fDigitsArray->IsSimulated()) {
1497 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1498 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1499 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1502 } // end of loop over time buckets
1503 } // end of lop over pads
1506 // This row has been digitized, delete nonused stuff
1509 for(lp=0;lp<nofDigits;lp++){
1510 if(pList[lp]) delete [] pList[lp];
1518 } // end of DigitizeRow
1520 //_____________________________________________________________________________
1522 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1523 TMatrixF *m1, TMatrixF *m2,Int_t *indexRange)
1526 //---------------------------------------------------------------
1527 // Calculates 2-D signal (pad,time) for a single track,
1528 // returns a pointer to the signal matrix and the track label
1529 // No digitization is performed at this level!!!
1530 //---------------------------------------------------------------
1532 //-----------------------------------------------------------------
1533 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1534 // Modified: Marian Ivanov
1535 //-----------------------------------------------------------------
1539 tv = (TVector*)p1->At(ntr); // pointer to a track
1542 Float_t label = v(0);
1543 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1)-1)/2;
1545 Int_t nElectrons = (tv->GetNrows()-1)/5;
1546 indexRange[0]=9999; // min pad
1547 indexRange[1]=-1; // max pad
1548 indexRange[2]=9999; //min time
1549 indexRange[3]=-1; // max time
1551 TMatrixF &signal = *m1;
1552 TMatrixF &total = *m2;
1554 // Loop over all electrons
1556 for(Int_t nel=0; nel<nElectrons; nel++){
1558 Float_t aval = v(idx+4);
1559 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1560 Float_t xyz[4]={v(idx+1),v(idx+2),v(idx+3),v(idx+5)};
1561 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
1563 Int_t *index = fTPCParam->GetResBin(0);
1564 Float_t *weight = & (fTPCParam->GetResWeight(0));
1566 if (n>0) for (Int_t i =0; i<n; i++){
1567 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1570 Int_t time=index[2];
1571 Float_t qweight = *(weight)*eltoadcfac;
1573 if (m1!=0) signal(pad,time)+=qweight;
1574 total(pad,time)+=qweight;
1575 if (indexRange[0]>pad) indexRange[0]=pad;
1576 if (indexRange[1]<pad) indexRange[1]=pad;
1577 if (indexRange[2]>time) indexRange[2]=time;
1578 if (indexRange[3]<time) indexRange[3]=time;
1585 } // end of loop over electrons
1587 return label; // returns track label when finished
1590 //_____________________________________________________________________________
1591 void AliTPC::GetList(Float_t label,Int_t np,TMatrixF *m,
1592 Int_t *indexRange, Float_t **pList)
1594 //----------------------------------------------------------------------
1595 // Updates the list of tracks contributing to digits for a given row
1596 //----------------------------------------------------------------------
1598 //-----------------------------------------------------------------
1599 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1600 //-----------------------------------------------------------------
1602 TMatrixF &signal = *m;
1604 // lop over nonzero digits
1606 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1607 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1610 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1612 if(signal(ip,it)<0.5) continue;
1614 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1616 if(!pList[globalIndex]){
1619 // Create new list (6 elements - 3 signals and 3 labels),
1622 pList[globalIndex] = new Float_t [6];
1626 *pList[globalIndex] = -1.;
1627 *(pList[globalIndex]+1) = -1.;
1628 *(pList[globalIndex]+2) = -1.;
1629 *(pList[globalIndex]+3) = -1.;
1630 *(pList[globalIndex]+4) = -1.;
1631 *(pList[globalIndex]+5) = -1.;
1633 *pList[globalIndex] = label;
1634 *(pList[globalIndex]+3) = signal(ip,it);
1638 // check the signal magnitude
1640 Float_t highest = *(pList[globalIndex]+3);
1641 Float_t middle = *(pList[globalIndex]+4);
1642 Float_t lowest = *(pList[globalIndex]+5);
1645 // compare the new signal with already existing list
1648 if(signal(ip,it)<lowest) continue; // neglect this track
1652 if (signal(ip,it)>highest){
1653 *(pList[globalIndex]+5) = middle;
1654 *(pList[globalIndex]+4) = highest;
1655 *(pList[globalIndex]+3) = signal(ip,it);
1657 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1658 *(pList[globalIndex]+1) = *pList[globalIndex];
1659 *pList[globalIndex] = label;
1661 else if (signal(ip,it)>middle){
1662 *(pList[globalIndex]+5) = middle;
1663 *(pList[globalIndex]+4) = signal(ip,it);
1665 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1666 *(pList[globalIndex]+1) = label;
1669 *(pList[globalIndex]+5) = signal(ip,it);
1670 *(pList[globalIndex]+2) = label;
1674 } // end of loop over pads
1675 } // end of loop over time bins
1678 //___________________________________________________________________
1679 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1680 Stat_t ntracks,TObjArray **row)
1683 //-----------------------------------------------------------------
1684 // Prepares the sector digitization, creates the vectors of
1685 // tracks for each row of this sector. The track vector
1686 // contains the track label and the position of electrons.
1687 //-----------------------------------------------------------------
1689 //-----------------------------------------------------------------
1690 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1691 //-----------------------------------------------------------------
1693 Float_t gasgain = fTPCParam->GetGasGain();
1697 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1700 if (fHitType>1) branch = TH->GetBranch("TPC2");
1701 else branch = TH->GetBranch("TPC");
1704 //----------------------------------------------
1705 // Create TObjArray-s, one for each row,
1706 // each TObjArray will store the TVectors
1707 // of electrons, one TVectors per each track.
1708 //----------------------------------------------
1710 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
1711 TVector **tracks = new TVector* [nrows+2]; //pointers to the track vectors
1713 for(i=0; i<nrows+2; i++){
1714 row[i] = new TObjArray;
1721 //--------------------------------------------------------------------
1722 // Loop over tracks, the "track" contains the full history
1723 //--------------------------------------------------------------------
1725 Int_t previousTrack,currentTrack;
1726 previousTrack = -1; // nothing to store so far!
1728 for(Int_t track=0;track<ntracks;track++){
1729 Bool_t isInSector=kTRUE;
1731 isInSector = TrackInVolume(isec,track);
1732 if (!isInSector) continue;
1734 branch->GetEntry(track); // get next track
1738 tpcHit = (AliTPChit*)FirstHit(-1);
1740 //--------------------------------------------------------------
1742 //--------------------------------------------------------------
1747 Int_t sector=tpcHit->fSector; // sector number
1749 tpcHit = (AliTPChit*) NextHit();
1753 currentTrack = tpcHit->Track(); // track number
1755 if(currentTrack != previousTrack){
1757 // store already filled fTrack
1759 for(i=0;i<nrows+2;i++){
1760 if(previousTrack != -1){
1761 if(nofElectrons[i]>0){
1762 TVector &v = *tracks[i];
1763 v(0) = previousTrack;
1764 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
1765 row[i]->Add(tracks[i]);
1768 delete tracks[i]; // delete empty TVector
1774 tracks[i] = new TVector(601); // TVectors for the next fTrack
1776 } // end of loop over rows
1778 previousTrack=currentTrack; // update track label
1781 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1783 //---------------------------------------------------
1784 // Calculate the electron attachment probability
1785 //---------------------------------------------------
1788 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
1789 /fTPCParam->GetDriftV();
1791 Float_t attProb = fTPCParam->GetAttCoef()*
1792 fTPCParam->GetOxyCont()*time; // fraction!
1794 //-----------------------------------------------
1795 // Loop over electrons
1796 //-----------------------------------------------
1799 for(Int_t nel=0;nel<qI;nel++){
1800 // skip if electron lost due to the attachment
1801 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1806 // protection for the nonphysical avalanche size (10**6 maximum)
1808 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
1809 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
1812 TransportElectron(xyz,index);
1814 fTPCParam->GetPadRow(xyz,index);
1815 // Electron track time (for pileup simulation)
1816 xyz[4] = tpcHit->Time()/fTPCParam->GetTSample();
1817 // row 0 - cross talk from the innermost row
1818 // row fNRow+1 cross talk from the outermost row
1819 rowNumber = index[2]+1;
1820 //transform position to local digit coordinates
1821 //relative to nearest pad row
1822 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
1824 if (isec <fTPCParam->GetNInnerSector()) {
1825 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
1826 y1 = fTPCParam->GetYInner(rowNumber);
1829 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
1830 y1 = fTPCParam->GetYOuter(rowNumber);
1832 // gain inefficiency at the wires edges - linear
1835 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));
1837 nofElectrons[rowNumber]++;
1838 //----------------------------------
1839 // Expand vector if necessary
1840 //----------------------------------
1841 if(nofElectrons[rowNumber]>120){
1842 Int_t range = tracks[rowNumber]->GetNrows();
1843 if((nofElectrons[rowNumber])>(range-1)/5){
1845 tracks[rowNumber]->ResizeTo(range+500); // Add 100 electrons
1849 TVector &v = *tracks[rowNumber];
1850 Int_t idx = 5*nofElectrons[rowNumber]-4;
1851 Real_t * position = &(((TVector&)v)(idx)); //make code faster
1852 memcpy(position,xyz,5*sizeof(Float_t));
1854 } // end of loop over electrons
1856 tpcHit = (AliTPChit*)NextHit();
1858 } // end of loop over hits
1859 } // end of loop over tracks
1862 // store remaining track (the last one) if not empty
1865 for(i=0;i<nrows+2;i++){
1866 if(nofElectrons[i]>0){
1867 TVector &v = *tracks[i];
1868 v(0) = previousTrack;
1869 tracks[i]->ResizeTo(5*nofElectrons[i]+1); // shrink if necessary
1870 row[i]->Add(tracks[i]);
1879 delete [] nofElectrons;
1881 } // end of MakeSector
1884 //_____________________________________________________________________________
1888 // Initialise TPC detector after definition of geometry
1890 AliDebug(1,"*********************************************");
1893 //_____________________________________________________________________________
1894 void AliTPC::MakeBranch(Option_t* option)
1897 // Create Tree branches for the TPC.
1900 Int_t buffersize = 4000;
1901 char branchname[10];
1902 sprintf(branchname,"%s",GetName());
1904 const char *h = strstr(option,"H");
1906 if ( h && (fHitType<=1) && (fHits == 0x0)) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
1908 AliDetector::MakeBranch(option);
1910 const char *d = strstr(option,"D");
1912 if (fDigits && fLoader->TreeD() && d) {
1913 MakeBranchInTree(gAlice->TreeD(), branchname, &fDigits, buffersize, 0);
1916 if (fHitType>1) MakeBranch2(option,0); // MI change 14.09.2000
1919 //_____________________________________________________________________________
1920 void AliTPC::ResetDigits()
1923 // Reset number of digits and the digits array for this detector
1926 if (fDigits) fDigits->Clear();
1929 //_____________________________________________________________________________
1930 void AliTPC::SetSecAL(Int_t sec)
1932 //---------------------------------------------------
1933 // Activate/deactivate selection for lower sectors
1934 //---------------------------------------------------
1936 //-----------------------------------------------------------------
1937 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1938 //-----------------------------------------------------------------
1942 //_____________________________________________________________________________
1943 void AliTPC::SetSecAU(Int_t sec)
1945 //----------------------------------------------------
1946 // Activate/deactivate selection for upper sectors
1947 //---------------------------------------------------
1949 //-----------------------------------------------------------------
1950 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1951 //-----------------------------------------------------------------
1955 //_____________________________________________________________________________
1956 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
1958 //----------------------------------------
1959 // Select active lower sectors
1960 //----------------------------------------
1962 //-----------------------------------------------------------------
1963 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1964 //-----------------------------------------------------------------
1974 //_____________________________________________________________________________
1975 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
1976 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
1977 Int_t s11 , Int_t s12)
1979 //--------------------------------
1980 // Select active upper sectors
1981 //--------------------------------
1983 //-----------------------------------------------------------------
1984 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1985 //-----------------------------------------------------------------
2001 //_____________________________________________________________________________
2002 void AliTPC::SetSens(Int_t sens)
2005 //-------------------------------------------------------------
2006 // Activates/deactivates the sensitive strips at the center of
2007 // the pad row -- this is for the space-point resolution calculations
2008 //-------------------------------------------------------------
2010 //-----------------------------------------------------------------
2011 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2012 //-----------------------------------------------------------------
2018 void AliTPC::SetSide(Float_t side=0.)
2020 // choice of the TPC side
2025 //____________________________________________________________________________
2026 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2027 Float_t p2,Float_t p3)
2030 // gax mixture definition
2042 //_____________________________________________________________________________
2044 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2047 // electron transport taking into account:
2049 // 2.ExB at the wires
2050 // 3. nonisochronity
2052 // xyz and index must be already transformed to system 1
2055 fTPCParam->Transform1to2(xyz,index);
2058 Float_t driftl=xyz[2];
2059 if(driftl<0.01) driftl=0.01;
2060 driftl=TMath::Sqrt(driftl);
2061 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2062 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2063 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2064 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2065 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2069 if (fTPCParam->GetMWPCReadout()==kTRUE){
2070 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2071 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2073 //add nonisochronity (not implemented yet)
2078 //_____________________________________________________________________________
2079 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2083 // Creates a TPC hit object
2094 //________________________________________________________________________
2095 // Additional code because of the AliTPCTrackHitsV2
2097 void AliTPC::MakeBranch2(Option_t *option,const char */*file*/)
2100 // Create a new branch in the current Root Tree
2101 // The branch of fHits is automatically split
2102 // MI change 14.09.2000
2104 if (fHitType<2) return;
2105 char branchname[10];
2106 sprintf(branchname,"%s2",GetName());
2108 // Get the pointer to the header
2109 const char *cH = strstr(option,"H");
2111 if (fTrackHits && TreeH() && cH && fHitType&4) {
2112 AliDebug(1,"Making branch for Type 4 Hits");
2113 TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2116 // if (fTrackHitsOld && TreeH() && cH && fHitType&2) {
2117 // AliDebug(1,"Making branch for Type 2 Hits");
2118 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2119 // TreeH(),fBufferSize,99);
2120 // TreeH()->GetListOfBranches()->Add(branch);
2124 void AliTPC::SetTreeAddress()
2126 //Sets tree address for hits
2128 if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2129 AliDetector::SetTreeAddress();
2131 if (fHitType>1) SetTreeAddress2();
2134 void AliTPC::SetTreeAddress2()
2137 // Set branch address for the TrackHits Tree
2142 char branchname[20];
2143 sprintf(branchname,"%s2",GetName());
2145 // Branch address for hit tree
2146 TTree *treeH = TreeH();
2147 if ((treeH)&&(fHitType&4)) {
2148 branch = treeH->GetBranch(branchname);
2150 branch->SetAddress(&fTrackHits);
2151 AliDebug(1,"fHitType&4 Setting");
2154 AliDebug(1,"fHitType&4 Failed (can not find branch)");
2157 // if ((treeH)&&(fHitType&2)) {
2158 // branch = treeH->GetBranch(branchname);
2160 // branch->SetAddress(&fTrackHitsOld);
2161 // AliDebug(1,"fHitType&2 Setting");
2164 // AliDebug(1,"fHitType&2 Failed (can not find branch)");
2166 //set address to TREETR
2168 TTree *treeTR = TreeTR();
2169 if (treeTR && fTrackReferences) {
2170 branch = treeTR->GetBranch(GetName());
2171 if (branch) branch->SetAddress(&fTrackReferences);
2176 void AliTPC::FinishPrimary()
2178 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2179 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2183 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2186 // add hit to the list
2189 int primary = gAlice->GetMCApp()->GetPrimary(track);
2190 gAlice->GetMCApp()->Particle(primary)->SetBit(kKeepBit);
2194 gAlice->GetMCApp()->FlagTrack(track);
2196 if (fTrackHits && fHitType&4)
2197 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2198 hits[1],hits[2],(Int_t)hits[3],hits[4]);
2199 // if (fTrackHitsOld &&fHitType&2 )
2200 // fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2201 // hits[1],hits[2],(Int_t)hits[3]);
2205 void AliTPC::ResetHits()
2207 if (fHitType&1) AliDetector::ResetHits();
2208 if (fHitType>1) ResetHits2();
2211 void AliTPC::ResetHits2()
2215 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2216 // if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2220 AliHit* AliTPC::FirstHit(Int_t track)
2222 if (fHitType>1) return FirstHit2(track);
2223 return AliDetector::FirstHit(track);
2225 AliHit* AliTPC::NextHit()
2230 if (fHitType>1) return NextHit2();
2232 return AliDetector::NextHit();
2235 AliHit* AliTPC::FirstHit2(Int_t track)
2238 // Initialise the hit iterator
2239 // Return the address of the first hit for track
2240 // If track>=0 the track is read from disk
2241 // while if track<0 the first hit of the current
2242 // track is returned
2245 gAlice->ResetHits();
2246 TreeH()->GetEvent(track);
2249 if (fTrackHits && fHitType&4) {
2250 fTrackHits->First();
2251 return fTrackHits->GetHit();
2253 // if (fTrackHitsOld && fHitType&2) {
2254 // fTrackHitsOld->First();
2255 // return fTrackHitsOld->GetHit();
2261 AliHit* AliTPC::NextHit2()
2264 //Return the next hit for the current track
2267 // if (fTrackHitsOld && fHitType&2) {
2268 // fTrackHitsOld->Next();
2269 // return fTrackHitsOld->GetHit();
2273 return fTrackHits->GetHit();
2279 void AliTPC::LoadPoints(Int_t)
2284 if(fHitType==1) AliDetector::LoadPoints(a);
2285 else LoadPoints2(a);
2289 void AliTPC::RemapTrackHitIDs(Int_t *map)
2294 if (!fTrackHits) return;
2296 // if (fTrackHitsOld && fHitType&2){
2297 // AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2298 // for (UInt_t i=0;i<arr->GetSize();i++){
2299 // AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2300 // info->fTrackID = map[info->fTrackID];
2303 // if (fTrackHitsOld && fHitType&4){
2304 if (fTrackHits && fHitType&4){
2305 TClonesArray * arr = fTrackHits->GetArray();;
2306 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2307 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2308 info->fTrackID = map[info->fTrackID];
2313 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2315 //return bool information - is track in given volume
2316 //load only part of the track information
2317 //return true if current track is in volume
2320 // if (fTrackHitsOld && fHitType&2) {
2321 // TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
2322 // br->GetEvent(track);
2323 // AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2324 // for (UInt_t j=0;j<ar->GetSize();j++){
2325 // if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2329 if (fTrackHits && fHitType&4) {
2330 TBranch * br1 = TreeH()->GetBranch("fVolumes");
2331 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
2332 br2->GetEvent(track);
2333 br1->GetEvent(track);
2334 Int_t *volumes = fTrackHits->GetVolumes();
2335 Int_t nvolumes = fTrackHits->GetNVolumes();
2336 if (!volumes && nvolumes>0) {
2337 AliWarning(Form("Problematic track\t%d\t%d",track,nvolumes));
2340 for (Int_t j=0;j<nvolumes; j++)
2341 if (volumes[j]==id) return kTRUE;
2345 TBranch * br = TreeH()->GetBranch("fSector");
2346 br->GetEvent(track);
2347 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2348 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2355 //_____________________________________________________________________________
2356 void AliTPC::LoadPoints2(Int_t)
2359 // Store x, y, z of all hits in memory
2361 // if (fTrackHits == 0 && fTrackHitsOld==0) return;
2362 if (fTrackHits == 0 ) return;
2365 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2366 // if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2368 if (nhits == 0) return;
2369 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2370 if (fPoints == 0) fPoints = new TObjArray(tracks);
2373 Int_t *ntrk=new Int_t[tracks];
2374 Int_t *limi=new Int_t[tracks];
2375 Float_t **coor=new Float_t*[tracks];
2376 for(Int_t i=0;i<tracks;i++) {
2382 AliPoints *points = 0;
2385 Int_t chunk=nhits/4+1;
2387 // Loop over all the hits and store their position
2389 ahit = FirstHit2(-1);
2391 trk=ahit->GetTrack();
2392 if(ntrk[trk]==limi[trk]) {
2394 // Initialise a new track
2395 fp=new Float_t[3*(limi[trk]+chunk)];
2397 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2398 delete [] coor[trk];
2405 fp[3*ntrk[trk] ] = ahit->X();
2406 fp[3*ntrk[trk]+1] = ahit->Y();
2407 fp[3*ntrk[trk]+2] = ahit->Z();
2415 for(trk=0; trk<tracks; ++trk) {
2417 points = new AliPoints();
2418 points->SetMarkerColor(GetMarkerColor());
2419 points->SetMarkerSize(GetMarkerSize());
2420 points->SetDetector(this);
2421 points->SetParticle(trk);
2422 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2423 fPoints->AddAt(points,trk);
2424 delete [] coor[trk];
2434 //_____________________________________________________________________________
2435 void AliTPC::LoadPoints3(Int_t)
2438 // Store x, y, z of all hits in memory
2439 // - only intersection point with pad row
2440 if (fTrackHits == 0) return;
2442 Int_t nhits = fTrackHits->GetEntriesFast();
2443 if (nhits == 0) return;
2444 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2445 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2446 fPoints->Expand(2*tracks);
2449 Int_t *ntrk=new Int_t[tracks];
2450 Int_t *limi=new Int_t[tracks];
2451 Float_t **coor=new Float_t*[tracks];
2452 for(Int_t i=0;i<tracks;i++) {
2458 AliPoints *points = 0;
2461 Int_t chunk=nhits/4+1;
2463 // Loop over all the hits and store their position
2465 ahit = FirstHit2(-1);
2469 trk=ahit->GetTrack();
2470 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2471 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2472 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2473 if (currentrow!=lastrow){
2474 lastrow = currentrow;
2475 //later calculate intersection point
2476 if(ntrk[trk]==limi[trk]) {
2478 // Initialise a new track
2479 fp=new Float_t[3*(limi[trk]+chunk)];
2481 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2482 delete [] coor[trk];
2489 fp[3*ntrk[trk] ] = ahit->X();
2490 fp[3*ntrk[trk]+1] = ahit->Y();
2491 fp[3*ntrk[trk]+2] = ahit->Z();
2498 for(trk=0; trk<tracks; ++trk) {
2500 points = new AliPoints();
2501 points->SetMarkerColor(GetMarkerColor()+1);
2502 points->SetMarkerStyle(5);
2503 points->SetMarkerSize(0.2);
2504 points->SetDetector(this);
2505 points->SetParticle(trk);
2506 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2507 fPoints->AddAt(points,tracks+trk);
2508 delete [] coor[trk];
2519 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
2522 fLoader = new AliTPCLoader(GetName(),topfoldername);
2526 ////////////////////////////////////////////////////////////////////////
2527 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
2529 // load TPC paarmeters from a given file or create new if the object
2530 // is not found there
2531 // 12/05/2003 This method should be moved to the AliTPCLoader
2532 // and one has to decide where to store the TPC parameters
2535 sprintf(paramName,"75x40_100x60_150x60");
2536 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
2538 AliDebugClass(1,Form("TPC parameters %s found.",paramName));
2540 AliWarningClass("TPC parameters not found. Create new (they may be incorrect)");
2541 paramTPC = new AliTPCParamSR;
2545 // the older version of parameters can be accessed with this code.
2546 // In some cases, we have old parameters saved in the file but
2547 // digits were created with new parameters, it can be distinguish
2548 // by the name of TPC TreeD. The code here is just for the case
2549 // we would need to compare with old data, uncomment it if needed.
2551 // char paramName[50];
2552 // sprintf(paramName,"75x40_100x60");
2553 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
2555 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2557 // sprintf(paramName,"75x40_100x60_150x60");
2558 // paramTPC=(AliTPCParam*)in->Get(paramName);
2560 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
2562 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
2564 // paramTPC = new AliTPCParamSR;