/**************************************************************************
* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
* *
* Author: The ALICE Off-line Project. *
* Contributors are mentioned in the code where appropriate. *
* *
* Permission to use, copy, modify and distribute this software and its *
* documentation strictly for non-commercial purposes is hereby granted *
* without fee, provided that the above copyright notice appears in all *
* copies and that both the copyright notice and this permission notice *
* appear in the supporting documentation. The authors make no claims *
* about the suitability of this software for any purpose. It is *
* provided "as is" without express or implied warranty. *
**************************************************************************/
/* $Id$ */
///////////////////////////////////////////////////////////////////////////////
// //
// Time Projection Chamber //
// This class contains the basic functions for the Time Projection Chamber //
// detector. Functions specific to one particular geometry are //
// contained in the derived classes //
// //
//Begin_Html
/*
*/
//End_Html
// //
// //
///////////////////////////////////////////////////////////////////////////////
//
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include "AliArrayBranch.h"
#include "AliClusters.h"
#include "AliComplexCluster.h"
#include "AliDigits.h"
#include "AliMagF.h"
#include "AliPoints.h"
#include "AliRun.h"
#include "AliRunLoader.h"
#include "AliSimDigits.h"
#include "AliTPC.h"
#include "AliTPC.h"
#include "AliTPCClustersArray.h"
#include "AliTPCClustersRow.h"
#include "AliTPCDigitsArray.h"
#include "AliTPCLoader.h"
#include "AliTPCPRF2D.h"
#include "AliTPCParamSR.h"
#include "AliTPCRF1D.h"
#include "AliTPCTrackHits.h"
#include "AliTPCTrackHitsV2.h"
#include "AliTPCcluster.h"
#include "AliTrackReference.h"
#include "AliMC.h"
#include "AliTPCDigitizer.h"
#include "AliTPCclustererMI.h"
#include "AliTPCtrackerMI.h"
#include "AliTPCpidESD.h"
ClassImp(AliTPC)
//_____________________________________________________________________________
// helper class for fast matrix and vector manipulation - no range checking
// origin - Marian Ivanov
class AliTPCFastMatrix : public TMatrix {
public :
AliTPCFastMatrix(Int_t rowlwb, Int_t rowupb, Int_t collwb, Int_t colupb);
#if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
Float_t & UncheckedAt(Int_t rown, Int_t coln) const {return fElements[(rown-fRowLwb)*fNcols+(coln-fColLwb)];} //fast acces
Float_t UncheckedAtFast(Int_t rown, Int_t coln) const {return fElements[(rown-fRowLwb)*fNcols+(coln-fColLwb)];} //fast acces
#else
Float_t & UncheckedAt(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
Float_t UncheckedAtFast(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
#endif
};
AliTPCFastMatrix::AliTPCFastMatrix(Int_t rowlwb, Int_t rowupb, Int_t collwb, Int_t colupb):
TMatrix(rowlwb, rowupb,collwb,colupb)
{
};
//_____________________________________________________________________________
class AliTPCFastVector : public TVector {
public :
AliTPCFastVector(Int_t size);
virtual ~AliTPCFastVector(){;}
Float_t & UncheckedAt(Int_t index) const {return fElements[index];} //fast acces
};
AliTPCFastVector::AliTPCFastVector(Int_t size):TVector(size){
};
//_____________________________________________________________________________
AliTPC::AliTPC()
{
//
// Default constructor
//
fIshunt = 0;
fHits = 0;
fDigits = 0;
fNsectors = 0;
fDigitsArray = 0;
fClustersArray = 0;
fDefaults = 0;
fTrackHits = 0;
fTrackHitsOld = 0;
fHitType = 2; //default CONTAINERS - based on ROOT structure
fTPCParam = 0;
fNoiseTable = 0;
fActiveSectors =0;
}
//_____________________________________________________________________________
AliTPC::AliTPC(const char *name, const char *title)
: AliDetector(name,title)
{
//
// Standard constructor
//
//
// Initialise arrays of hits and digits
fHits = new TClonesArray("AliTPChit", 176);
gAlice->GetMCApp()->AddHitList(fHits);
fDigitsArray = 0;
fClustersArray= 0;
fDefaults = 0;
//
fTrackHits = new AliTPCTrackHitsV2;
fTrackHits->SetHitPrecision(0.002);
fTrackHits->SetStepPrecision(0.003);
fTrackHits->SetMaxDistance(100);
fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
fTrackHitsOld->SetHitPrecision(0.002);
fTrackHitsOld->SetStepPrecision(0.003);
fTrackHitsOld->SetMaxDistance(100);
fNoiseTable =0;
fHitType = 2;
fActiveSectors = 0;
//
// Initialise counters
fNsectors = 0;
//
fIshunt = 0;
//
// Initialise color attributes
SetMarkerColor(kYellow);
//
// Set TPC parameters
//
if (!strcmp(title,"Default")) {
fTPCParam = new AliTPCParamSR;
} else {
cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
fTPCParam=0;
}
}
//_____________________________________________________________________________
AliTPC::AliTPC(const AliTPC& t):AliDetector(t){
//
// dummy copy constructor
//
}
AliTPC::~AliTPC()
{
//
// TPC destructor
//
fIshunt = 0;
delete fHits;
delete fDigits;
delete fTPCParam;
delete fTrackHits; //MI 15.09.2000
delete fTrackHitsOld; //MI 10.12.2001
if (fNoiseTable) delete [] fNoiseTable;
}
//_____________________________________________________________________________
void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
{
//
// Add a hit to the list
//
// TClonesArray &lhits = *fHits;
// new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
if (fHitType&1){
TClonesArray &lhits = *fHits;
new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
}
if (fHitType>1)
AddHit2(track,vol,hits);
}
//_____________________________________________________________________________
void AliTPC::BuildGeometry()
{
//
// Build TPC ROOT TNode geometry for the event display
//
TNode *nNode, *nTop;
TTUBS *tubs;
Int_t i;
const int kColorTPC=19;
char name[5], title[25];
const Double_t kDegrad=TMath::Pi()/180;
const Double_t kRaddeg=180./TMath::Pi();
Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
Int_t nLo = fTPCParam->GetNInnerSector()/2;
Int_t nHi = fTPCParam->GetNOuterSector()/2;
const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
Double_t rl,ru;
//
// Get ALICE top node
//
nTop=gAlice->GetGeometry()->GetNode("alice");
// inner sectors
rl = fTPCParam->GetInnerRadiusLow();
ru = fTPCParam->GetInnerRadiusUp();
for(i=0;iSetNumberOfDivisions(1);
nTop->cd();
nNode = new TNode(name,title,name,0,0,0,"");
nNode->SetLineColor(kColorTPC);
fNodes->Add(nNode);
}
// Outer sectors
rl = fTPCParam->GetOuterRadiusLow();
ru = fTPCParam->GetOuterRadiusUp();
for(i=0;iSetNumberOfDivisions(1);
nTop->cd();
nNode = new TNode(name,title,name,0,0,0,"");
nNode->SetLineColor(kColorTPC);
fNodes->Add(nNode);
}
}
//_____________________________________________________________________________
Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t ) const
{
//
// Calculate distance from TPC to mouse on the display
// Dummy procedure
//
return 9999;
}
void AliTPC::Clusters2Tracks() const
{
//-----------------------------------------------------------------
// This is a track finder.
//-----------------------------------------------------------------
Error("Clusters2Tracks",
"Dummy function ! Call AliTPCtracker::Clusters2Tracks(...) instead !");
}
//_____________________________________________________________________________
void AliTPC::Reconstruct() const
{
// reconstruct clusters
AliLoader* loader = GetLoader();
loader->LoadRecPoints("recreate");
loader->LoadDigits("read");
AliTPCclustererMI clusterer(fTPCParam);
AliRunLoader* runLoader = loader->GetRunLoader();
Int_t nEvents = runLoader->GetNumberOfEvents();
for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
runLoader->GetEvent(iEvent);
TTree* treeClusters = loader->TreeR();
if (!treeClusters) {
loader->MakeTree("R");
treeClusters = loader->TreeR();
}
TTree* treeDigits = loader->TreeD();
if (!treeDigits) {
Error("Reconstruct", "Can't get digits tree !");
return;
}
// clusterer.Digits2Clusters(treeDigits, treeClusters);
clusterer.SetInput(treeDigits);
clusterer.SetOutput(treeClusters);
clusterer.Digits2Clusters();
loader->WriteRecPoints("OVERWRITE");
}
loader->UnloadRecPoints();
loader->UnloadDigits();
}
//_____________________________________________________________________________
AliTracker* AliTPC::CreateTracker() const
{
// create a TPC tracker
return new AliTPCtrackerMI(fTPCParam);
}
//_____________________________________________________________________________
void AliTPC::FillESD(AliESD* esd) const
{
// make PID
Double_t parTPC[] = {47., 0.10, 10.};
AliTPCpidESD tpcPID(parTPC);
tpcPID.MakePID(esd);
}
//_____________________________________________________________________________
void AliTPC::CreateMaterials()
{
//-----------------------------------------------
// Create Materials for for TPC simulations
//-----------------------------------------------
//-----------------------------------------------------------------
// Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
//-----------------------------------------------------------------
Int_t iSXFLD=gAlice->Field()->Integ();
Float_t sXMGMX=gAlice->Field()->Max();
Float_t amat[5]; // atomic numbers
Float_t zmat[5]; // z
Float_t wmat[5]; // proportions
Float_t density;
Float_t apure[2];
//***************** Gases *************************
//-------------------------------------------------
// pure gases
//-------------------------------------------------
// Neon
amat[0]= 20.18;
zmat[0]= 10.;
density = 0.0009;
apure[0]=amat[0];
AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
// Argon
amat[0]= 39.948;
zmat[0]= 18.;
density = 0.001782;
apure[1]=amat[0];
AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
//--------------------------------------------------------------
// gases - compounds
//--------------------------------------------------------------
Float_t amol[3];
// CO2
amat[0]=12.011;
amat[1]=15.9994;
zmat[0]=6.;
zmat[1]=8.;
wmat[0]=1.;
wmat[1]=2.;
density=0.001977;
amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
// CF4
amat[0]=12.011;
amat[1]=18.998;
zmat[0]=6.;
zmat[1]=9.;
wmat[0]=1.;
wmat[1]=4.;
density=0.003034;
amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
// CH4
amat[0]=12.011;
amat[1]=1.;
zmat[0]=6.;
zmat[1]=1.;
wmat[0]=1.;
wmat[1]=4.;
density=0.000717;
amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
//----------------------------------------------------------------
// gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
//----------------------------------------------------------------
char namate[21]="";
density = 0.;
Float_t am=0;
Int_t nc;
Float_t rho,absl,x0,buf[1];
Int_t nbuf;
Float_t a,z;
for(nc = 0;ncGfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,x0,absl,buf,nbuf);
amat[nc] = a;
zmat[nc] = z;
Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
density += fMixtProp[nc]*rho; // density of the mixture
}
// mixture proportions by weight!
for(nc = 0;nc=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
apure[nnc] : amol[nnc])/am;
}
// Drift gases 1 - nonsensitive, 2 - sensitive
AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
// Air
amat[0] = 14.61;
zmat[0] = 7.3;
density = 0.001205;
AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
//----------------------------------------------------------------------
// solid materials
//----------------------------------------------------------------------
// Kevlar C14H22O2N2
amat[0] = 12.011;
amat[1] = 1.;
amat[2] = 15.999;
amat[3] = 14.006;
zmat[0] = 6.;
zmat[1] = 1.;
zmat[2] = 8.;
zmat[3] = 7.;
wmat[0] = 14.;
wmat[1] = 22.;
wmat[2] = 2.;
wmat[3] = 2.;
density = 1.45;
AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
// NOMEX
amat[0] = 12.011;
amat[1] = 1.;
amat[2] = 15.999;
amat[3] = 14.006;
zmat[0] = 6.;
zmat[1] = 1.;
zmat[2] = 8.;
zmat[3] = 7.;
wmat[0] = 14.;
wmat[1] = 22.;
wmat[2] = 2.;
wmat[3] = 2.;
density = 0.03;
AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
// Makrolon C16H18O3
amat[0] = 12.011;
amat[1] = 1.;
amat[2] = 15.999;
zmat[0] = 6.;
zmat[1] = 1.;
zmat[2] = 8.;
wmat[0] = 16.;
wmat[1] = 18.;
wmat[2] = 3.;
density = 1.2;
AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
// Mylar C5H4O2
amat[0]=12.011;
amat[1]=1.;
amat[2]=15.9994;
zmat[0]=6.;
zmat[1]=1.;
zmat[2]=8.;
wmat[0]=5.;
wmat[1]=4.;
wmat[2]=2.;
density = 1.39;
AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
// SiO2 - used later for the glass fiber
amat[0]=28.086;
amat[1]=15.9994;
zmat[0]=14.;
zmat[1]=8.;
wmat[0]=1.;
wmat[1]=2.;
AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
// Al
amat[0] = 26.98;
zmat[0] = 13.;
density = 2.7;
AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
// Si
amat[0] = 28.086;
zmat[0] = 14.;
density = 2.33;
AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
// Cu
amat[0] = 63.546;
zmat[0] = 29.;
density = 8.96;
AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
// Tedlar C2H3F
amat[0] = 12.011;
amat[1] = 1.;
amat[2] = 18.998;
zmat[0] = 6.;
zmat[1] = 1.;
zmat[2] = 9.;
wmat[0] = 2.;
wmat[1] = 3.;
wmat[2] = 1.;
density = 1.71;
AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
// Plexiglas C5H8O2
amat[0]=12.011;
amat[1]=1.;
amat[2]=15.9994;
zmat[0]=6.;
zmat[1]=1.;
zmat[2]=8.;
wmat[0]=5.;
wmat[1]=8.;
wmat[2]=2.;
density=1.18;
AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
// Epoxy - C14 H20 O3
amat[0]=12.011;
amat[1]=1.;
amat[2]=15.9994;
zmat[0]=6.;
zmat[1]=1.;
zmat[2]=8.;
wmat[0]=14.;
wmat[1]=20.;
wmat[2]=3.;
density=1.25;
AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
// Carbon
amat[0]=12.011;
zmat[0]=6.;
density= 2.265;
AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
// get epoxy
gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,x0,absl,buf,nbuf);
// Carbon fiber
wmat[0]=0.644; // by weight!
wmat[1]=0.356;
density=0.5*(1.25+2.265);
AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
// get SiO2
gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,x0,absl,buf,nbuf);
wmat[0]=0.725; // by weight!
wmat[1]=0.275;
density=1.7;
AliMixture(39,"G10",amat,zmat,density,2,wmat);
//----------------------------------------------------------
// tracking media for gases
//----------------------------------------------------------
AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
//-----------------------------------------------------------
// tracking media for solids
//-----------------------------------------------------------
AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
}
void AliTPC::GenerNoise(Int_t tablesize)
{
//
//Generate table with noise
//
if (fTPCParam==0) {
// error message
fNoiseDepth=0;
return;
}
if (fNoiseTable) delete[] fNoiseTable;
fNoiseTable = new Float_t[tablesize];
fNoiseDepth = tablesize;
fCurrentNoise =0; //!index of the noise in the noise table
Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
for (Int_t i=0;iGaus(0,norm);
}
Float_t AliTPC::GetNoise()
{
// get noise from table
// if ((fCurrentNoise%10)==0)
// fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
return fNoiseTable[fCurrentNoise++];
//gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
}
Bool_t AliTPC::IsSectorActive(Int_t sec) const
{
//
// check if the sector is active
if (!fActiveSectors) return kTRUE;
else return fActiveSectors[sec];
}
void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
{
// activate interesting sectors
SetTreeAddress();//just for security
if (fActiveSectors) delete [] fActiveSectors;
fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
for (Int_t i=0;iGetNSector();i++) fActiveSectors[i]=kFALSE;
for (Int_t i=0;i=0) && sectors[i]GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
}
void AliTPC::SetActiveSectors(Int_t flag)
{
//
// activate sectors which were hitted by tracks
//loop over tracks
SetTreeAddress();//just for security
if (fHitType==0) return; // if Clones hit - not short volume ID information
if (fActiveSectors) delete [] fActiveSectors;
fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
if (flag) {
for (Int_t i=0;iGetNSector();i++) fActiveSectors[i]=kTRUE;
return;
}
for (Int_t i=0;iGetNSector();i++) fActiveSectors[i]=kFALSE;
TBranch * branch=0;
if (TreeH() == 0x0)
{
Fatal("SetActiveSectors","Can not find TreeH in folder");
return;
}
if (fHitType>1) branch = TreeH()->GetBranch("TPC2");
else branch = TreeH()->GetBranch("TPC");
Stat_t ntracks = TreeH()->GetEntries();
// loop over all hits
if (GetDebug()) cout<<"\nAliTPC::SetActiveSectors(): Got "<GetBranch("fVolumes");
TBranch * br2 = TreeH()->GetBranch("fNVolumes");
br1->GetEvent(track);
br2->GetEvent(track);
Int_t *volumes = fTrackHits->GetVolumes();
for (Int_t j=0;jGetNVolumes(); j++)
fActiveSectors[volumes[j]]=kTRUE;
}
//
if (fTrackHitsOld && fHitType&2) {
TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
br->GetEvent(track);
AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
for (UInt_t j=0;jGetSize();j++){
fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
}
}
}
}
void AliTPC::Digits2Clusters(Int_t /*eventnumber*/) const
{
//-----------------------------------------------------------------
// This is a simple cluster finder.
//-----------------------------------------------------------------
Error("Digits2Clusters",
"Dummy function ! Call AliTPCclusterer::Digits2Clusters(...) instead !");
}
extern Double_t SigmaY2(Double_t, Double_t, Double_t);
extern Double_t SigmaZ2(Double_t, Double_t);
//_____________________________________________________________________________
void AliTPC::Hits2Clusters(Int_t /*eventn*/)
{
//--------------------------------------------------------
// TPC simple cluster generator from hits
// obtained from the TPC Fast Simulator
// The point errors are taken from the parametrization
//--------------------------------------------------------
//-----------------------------------------------------------------
// Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
//-----------------------------------------------------------------
// Adopted to Marian's cluster data structure by I.Belikov, CERN,
// Jouri.Belikov@cern.ch
//----------------------------------------------------------------
/////////////////////////////////////////////////////////////////////////////
//
//---------------------------------------------------------------------
// ALICE TPC Cluster Parameters
//--------------------------------------------------------------------
// Cluster width in rphi
const Float_t kACrphi=0.18322;
const Float_t kBCrphi=0.59551e-3;
const Float_t kCCrphi=0.60952e-1;
// Cluster width in z
const Float_t kACz=0.19081;
const Float_t kBCz=0.55938e-3;
const Float_t kCCz=0.30428;
if (!fLoader) {
cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
return;
}
//if(fDefaults == 0) SetDefaults();
Float_t sigmaRphi,sigmaZ,clRphi,clZ;
//
TParticle *particle; // pointer to a given particle
AliTPChit *tpcHit; // pointer to a sigle TPC hit
Int_t sector;
Int_t ipart;
Float_t xyz[5];
Float_t pl,pt,tanth,rpad,ratio;
Float_t cph,sph;
//---------------------------------------------------------------
// Get the access to the tracks
//---------------------------------------------------------------
TTree *tH = TreeH();
if (tH == 0x0)
{
Fatal("Hits2Clusters","Can not find TreeH in folder");
return;
}
SetTreeAddress();
Stat_t ntracks = tH->GetEntries();
//Switch to the output file
if (fLoader->TreeR() == 0x0) fLoader->MakeTree("R");
cout<<"fTPCParam->GetTitle() = "<GetTitle()<GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
rl->CdGAFile();
//fTPCParam->Write(fTPCParam->GetTitle());
AliTPCClustersArray carray;
carray.Setup(fTPCParam);
carray.SetClusterType("AliTPCcluster");
carray.MakeTree(fLoader->TreeR());
Int_t nclusters=0; //cluster counter
//------------------------------------------------------------
// Loop over all sectors (72 sectors for 20 deg
// segmentation for both lower and upper sectors)
// Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
// Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
//
// First cluster for sector 0 starts at "0"
//------------------------------------------------------------
for(Int_t isec=0;isecGetNSector();isec++){
//MI change
fTPCParam->AdjustCosSin(isec,cph,sph);
//------------------------------------------------------------
// Loop over tracks
//------------------------------------------------------------
for(Int_t track=0;trackGetEvent(track);
//
// Get number of the TPC hits
//
tpcHit = (AliTPChit*)FirstHit(-1);
// Loop over hits
//
while(tpcHit){
if (tpcHit->fQ == 0.) {
tpcHit = (AliTPChit*) NextHit();
continue; //information about track (I.Belikov)
}
sector=tpcHit->fSector; // sector number
if(sector != isec){
tpcHit = (AliTPChit*) NextHit();
continue;
}
ipart=tpcHit->Track();
particle=gAlice->GetMCApp()->Particle(ipart);
pl=particle->Pz();
pt=particle->Pt();
if(pt < 1.e-9) pt=1.e-9;
tanth=pl/pt;
tanth = TMath::Abs(tanth);
rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
// space-point resolutions
sigmaRphi=SigmaY2(rpad,tanth,pt);
sigmaZ =SigmaZ2(rpad,tanth );
// cluster widths
clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
// temporary protection
if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
if(sigmaZ < 0.) sigmaZ=0.4e-3;
if(clRphi < 0.) clRphi=2.5e-3;
if(clZ < 0.) clZ=2.5e-5;
//
//
// smearing --> rotate to the 1 (13) or to the 25 (49) sector,
// then the inaccuracy in a X-Y plane is only along Y (pad row)!
//
Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
Float_t ymax=xprim*TMath::Tan(0.5*alpha);
if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
xyz[2]=sigmaRphi; // fSigmaY2
xyz[3]=sigmaZ; // fSigmaZ2
xyz[4]=tpcHit->fQ; // q
AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
Int_t tracks[3]={tpcHit->Track(), -1, -1};
AliTPCcluster cluster(tracks,xyz);
clrow->InsertCluster(&cluster); nclusters++;
tpcHit = (AliTPChit*)NextHit();
} // end of loop over hits
} // end of loop over tracks
Int_t nrows=fTPCParam->GetNRow(isec);
for (Int_t irow=0; irowWriteRecPoints("OVERWRITE");
} // end of function
//_________________________________________________________________
void AliTPC::Hits2ExactClustersSector(Int_t isec)
{
//--------------------------------------------------------
//calculate exact cross point of track and given pad row
//resulting values are expressed in "digit" coordinata
//--------------------------------------------------------
//-----------------------------------------------------------------
// Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
//-----------------------------------------------------------------
//
if (fClustersArray==0){
return;
}
//
TParticle *particle; // pointer to a given particle
AliTPChit *tpcHit; // pointer to a sigle TPC hit
// Int_t sector,nhits;
Int_t ipart;
const Int_t kcmaxhits=30000;
AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
AliTPCFastVector & xxx = *xxxx;
Int_t maxhits = kcmaxhits;
//construct array for each padrow
for (Int_t i=0; iGetNRow(isec);i++)
fClustersArray->CreateRow(isec,i);
//---------------------------------------------------------------
// Get the access to the tracks
//---------------------------------------------------------------
TTree *tH = TreeH();
if (tH == 0x0)
{
Fatal("Hits2Clusters","Can not find TreeH in folder");
return;
}
SetTreeAddress();
Stat_t ntracks = tH->GetEntries();
Int_t npart = gAlice->GetMCApp()->GetNtrack();
//MI change
TBranch * branch=0;
if (fHitType>1) branch = tH->GetBranch("TPC2");
else branch = tH->GetBranch("TPC");
//------------------------------------------------------------
// Loop over tracks
//------------------------------------------------------------
for(Int_t track=0;trackGetEntry(track); // get next track
//
// Get number of the TPC hits and a pointer
// to the particles
// Loop over hits
//
Int_t currentIndex=0;
Int_t lastrow=-1; //last writen row
//M.I. changes
tpcHit = (AliTPChit*)FirstHit(-1);
while(tpcHit){
Int_t sector=tpcHit->fSector; // sector number
if(sector != isec){
tpcHit = (AliTPChit*) NextHit();
continue;
}
ipart=tpcHit->Track();
if (ipartGetMCApp()->Particle(ipart);
//find row number
Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
Int_t index[3]={1,isec,0};
Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
if (lastrow<0) lastrow=currentrow;
if (currentrow==lastrow){
if ( currentIndex>=maxhits){
maxhits+=kcmaxhits;
xxx.ResizeTo(4*maxhits);
}
xxx(currentIndex*4)=x[0];
xxx(currentIndex*4+1)=x[1];
xxx(currentIndex*4+2)=x[2];
xxx(currentIndex*4+3)=tpcHit->fQ;
currentIndex++;
}
else
if (currentIndex>2){
Float_t sumx=0;
Float_t sumx2=0;
Float_t sumx3=0;
Float_t sumx4=0;
Float_t sumy=0;
Float_t sumxy=0;
Float_t sumx2y=0;
Float_t sumz=0;
Float_t sumxz=0;
Float_t sumx2z=0;
Float_t sumq=0;
for (Int_t index=0;indexGetNPads(isec,lastrow)-1)/2;
Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
sumx2*(sumx*sumx3-sumx2*sumx2);
Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
sumx2*(sumxy*sumx3-sumx2y*sumx2);
Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
sumx2*(sumxz*sumx3-sumx2z*sumx2);
Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
sumx2*(sumx*sumx2y-sumx2*sumxy);
Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
sumx2*(sumx*sumx2z-sumx2*sumxz);
if (TMath::Abs(det)<0.00001){
tpcHit = (AliTPChit*)NextHit();
continue;
}
Float_t y=detay/det+centralPad;
Float_t z=detaz/det;
Float_t by=detby/det; //y angle
Float_t bz=detbz/det; //z angle
sumy/=Float_t(currentIndex);
sumz/=Float_t(currentIndex);
AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
if (row!=0) {
AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
// AliComplexCluster cl;
cl->fX=z;
cl->fY=y;
cl->fQ=sumq;
cl->fSigmaX2=bz;
cl->fSigmaY2=by;
cl->fTracks[0]=ipart;
}
currentIndex=0;
lastrow=currentrow;
} //end of calculating cluster for given row
tpcHit = (AliTPChit*)NextHit();
} // end of loop over hits
} // end of loop over tracks
//write padrows to tree
for (Int_t ii=0; iiGetNRow(isec);ii++) {
fClustersArray->StoreRow(isec,ii);
fClustersArray->ClearRow(isec,ii);
}
xxxx->Delete();
}
//______________________________________________________________________
AliDigitizer* AliTPC::CreateDigitizer(AliRunDigitizer* manager) const
{
return new AliTPCDigitizer(manager);
}
//__
void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
{
//create digits from summable digits
GenerNoise(500000); //create teble with noise
//conect tree with sSDigits
TTree *t = fLoader->TreeS();
if (t == 0x0)
{
fLoader->LoadSDigits("READ");
t = fLoader->TreeS();
if (t == 0x0)
{
Error("SDigits2Digits2","Can not get input TreeS");
return;
}
}
if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
AliSimDigits digarr, *dummy=&digarr;
TBranch* sdb = t->GetBranch("Segment");
if (sdb == 0x0)
{
Error("SDigits2Digits2","Can not find branch with segments in TreeS.");
return;
}
sdb->SetAddress(&dummy);
Stat_t nentries = t->GetEntries();
// set zero suppression
fTPCParam->SetZeroSup(2);
// get zero suppression
Int_t zerosup = fTPCParam->GetZeroSup();
//make tree with digits
AliTPCDigitsArray *arr = new AliTPCDigitsArray;
arr->SetClass("AliSimDigits");
arr->Setup(fTPCParam);
arr->MakeTree(fLoader->TreeD());
AliTPCParam * par = fTPCParam;
//Loop over segments of the TPC
for (Int_t n=0; nGetEvent(n);
Int_t sec, row;
if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
cerr<<"AliTPC warning: invalid segment ID ! "<CreateRow(sec,row);
Int_t nrows = digrow->GetNRows();
Int_t ncols = digrow->GetNCols();
digrow->ExpandBuffer();
digarr.ExpandBuffer();
digrow->ExpandTrackBuffer();
digarr.ExpandTrackBuffer();
Short_t * pamp0 = digarr.GetDigits();
Int_t * ptracks0 = digarr.GetTracks();
Short_t * pamp1 = digrow->GetDigits();
Int_t * ptracks1 = digrow->GetTracks();
Int_t nelems =nrows*ncols;
Int_t saturation = fTPCParam->GetADCSat();
//use internal structure of the AliDigits - for speed reason
//if you cahnge implementation
//of the Alidigits - it must be rewriten -
for (Int_t i= 0; izerosup){
if (q>saturation) q=saturation;
*pamp1=(Short_t)q;
//if (ptracks0[0]==0)
// ptracks1[0]=1;
//else
ptracks1[0]=ptracks0[0];
ptracks1[nelems]=ptracks0[nelems];
ptracks1[2*nelems]=ptracks0[2*nelems];
}
pamp0++;
pamp1++;
ptracks0++;
ptracks1++;
}
arr->StoreRow(sec,row);
arr->ClearRow(sec,row);
// cerr<WriteDigits("OVERWRITE");
delete arr;
}
//__________________________________________________________________
void AliTPC::SetDefaults(){
//
// setting the defaults
//
// cerr<<"Setting default parameters...\n";
// Set response functions
//
AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
rl->CdGAFile();
AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
if(param){
printf("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...\n");
delete param;
param = new AliTPCParamSR();
}
else {
param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
}
if(!param){
printf("No TPC parameters found\n");
exit(4);
}
AliTPCPRF2D * prfinner = new AliTPCPRF2D;
AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
rf->SetOffset(3*param->GetZSigma());
rf->Update();
TDirectory *savedir=gDirectory;
TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
if (!f->IsOpen()) {
cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
exit(3);
}
TString s;
prfinner->Read("prf_07504_Gati_056068_d02");
//PH Set different names
s=prfinner->GetGRF()->GetName();
s+="in";
prfinner->GetGRF()->SetName(s.Data());
prfouter1->Read("prf_10006_Gati_047051_d03");
s=prfouter1->GetGRF()->GetName();
s+="out1";
prfouter1->GetGRF()->SetName(s.Data());
prfouter2->Read("prf_15006_Gati_047051_d03");
s=prfouter2->GetGRF()->GetName();
s+="out2";
prfouter2->GetGRF()->SetName(s.Data());
f->Close();
savedir->cd();
param->SetInnerPRF(prfinner);
param->SetOuter1PRF(prfouter1);
param->SetOuter2PRF(prfouter2);
param->SetTimeRF(rf);
// set fTPCParam
SetParam(param);
fDefaults = 1;
}
//__________________________________________________________________
void AliTPC::Hits2Digits()
{
//
// creates digits from hits
//
fLoader->LoadHits("read");
fLoader->LoadDigits("recreate");
AliRunLoader* runLoader = fLoader->GetRunLoader();
for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
runLoader->GetEvent(iEvent);
SetActiveSectors();
Hits2Digits(iEvent);
}
fLoader->UnloadHits();
fLoader->UnloadDigits();
}
//__________________________________________________________________
void AliTPC::Hits2Digits(Int_t eventnumber)
{
//----------------------------------------------------
// Loop over all sectors for a single event
//----------------------------------------------------
AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
rl->GetEvent(eventnumber);
if (fLoader->TreeH() == 0x0)
{
if(fLoader->LoadHits())
{
Error("Hits2Digits","Can not load hits.");
}
}
SetTreeAddress();
if (fLoader->TreeD() == 0x0 )
{
fLoader->MakeTree("D");
if (fLoader->TreeD() == 0x0 )
{
Error("Hits2Digits","Can not get TreeD");
return;
}
}
if(fDefaults == 0) SetDefaults(); // check if the parameters are set
GenerNoise(500000); //create teble with noise
//setup TPCDigitsArray
if(GetDigitsArray()) delete GetDigitsArray();
AliTPCDigitsArray *arr = new AliTPCDigitsArray;
arr->SetClass("AliSimDigits");
arr->Setup(fTPCParam);
arr->MakeTree(fLoader->TreeD());
SetDigitsArray(arr);
fDigitsSwitch=0; // standard digits
// cerr<<"Digitizing TPC -- normal digits...\n";
for(Int_t isec=0;isecGetNSector();isec++)
if (IsSectorActive(isec))
{
if (fDebug) Info("Hits2Digits","Sector %d is active.",isec);
Hits2DigitsSector(isec);
}
else
{
if (fDebug) Info("Hits2Digits","Sector %d is NOT active.",isec);
}
fLoader->WriteDigits("OVERWRITE");
//this line prevents the crash in the similar one
//on the beginning of this method
//destructor attempts to reset the tree, which is deleted by the loader
//need to be redesign
if(GetDigitsArray()) delete GetDigitsArray();
SetDigitsArray(0x0);
}
//__________________________________________________________________
void AliTPC::Hits2SDigits2(Int_t eventnumber)
{
//-----------------------------------------------------------
// summable digits - 16 bit "ADC", no noise, no saturation
//-----------------------------------------------------------
//----------------------------------------------------
// Loop over all sectors for a single event
//----------------------------------------------------
// AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
AliRunLoader* rl = fLoader->GetRunLoader();
rl->GetEvent(eventnumber);
if (fLoader->TreeH() == 0x0)
{
if(fLoader->LoadHits())
{
Error("Hits2Digits","Can not load hits.");
return;
}
}
SetTreeAddress();
if (fLoader->TreeS() == 0x0 )
{
fLoader->MakeTree("S");
}
if(fDefaults == 0) SetDefaults();
GenerNoise(500000); //create table with noise
//setup TPCDigitsArray
if(GetDigitsArray()) delete GetDigitsArray();
AliTPCDigitsArray *arr = new AliTPCDigitsArray;
arr->SetClass("AliSimDigits");
arr->Setup(fTPCParam);
arr->MakeTree(fLoader->TreeS());
SetDigitsArray(arr);
// cerr<<"Digitizing TPC -- summable digits...\n";
fDigitsSwitch=1; // summable digits
// set zero suppression to "0"
fTPCParam->SetZeroSup(0);
for(Int_t isec=0;isecGetNSector();isec++)
if (IsSectorActive(isec))
{
// cout<<"Sector "<WriteSDigits("OVERWRITE");
//this line prevents the crash in the similar one
//on the beginning of this method
//destructor attempts to reset the tree, which is deleted by the loader
//need to be redesign
if(GetDigitsArray()) delete GetDigitsArray();
SetDigitsArray(0x0);
}
//__________________________________________________________________
void AliTPC::Hits2SDigits()
{
//-----------------------------------------------------------
// summable digits - 16 bit "ADC", no noise, no saturation
//-----------------------------------------------------------
fLoader->LoadHits("read");
fLoader->LoadSDigits("recreate");
AliRunLoader* runLoader = fLoader->GetRunLoader();
for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
runLoader->GetEvent(iEvent);
SetTreeAddress();
SetActiveSectors();
Hits2SDigits2(iEvent);
}
fLoader->UnloadHits();
fLoader->UnloadSDigits();
}
//_____________________________________________________________________________
void AliTPC::Hits2DigitsSector(Int_t isec)
{
//-------------------------------------------------------------------
// TPC conversion from hits to digits.
//-------------------------------------------------------------------
//-----------------------------------------------------------------
// Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
//-----------------------------------------------------------------
//-------------------------------------------------------
// Get the access to the track hits
//-------------------------------------------------------
// check if the parameters are set - important if one calls this method
// directly, not from the Hits2Digits
if(fDefaults == 0) SetDefaults();
TTree *tH = TreeH(); // pointer to the hits tree
if (tH == 0x0)
{
Fatal("Hits2DigitsSector","Can not find TreeH in folder");
return;
}
Stat_t ntracks = tH->GetEntries();
if( ntracks > 0){
//-------------------------------------------
// Only if there are any tracks...
//-------------------------------------------
TObjArray **row;
//printf("*** Processing sector number %d ***\n",isec);
Int_t nrows =fTPCParam->GetNRow(isec);
row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
MakeSector(isec,nrows,tH,ntracks,row);
//--------------------------------------------------------
// Digitize this sector, row by row
// row[i] is the pointer to the TObjArray of AliTPCFastVectors,
// each one containing electrons accepted on this
// row, assigned into tracks
//--------------------------------------------------------
Int_t i;
if (fDigitsArray->GetTree()==0)
{
Fatal("Hits2DigitsSector","Tree not set in fDigitsArray");
}
for (i=0;iCreateRow(isec,i);
DigitizeRow(i,isec,row);
fDigitsArray->StoreRow(isec,i);
Int_t ndig = dig->GetDigitSize();
if (gDebug > 10)
printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
fDigitsArray->ClearRow(isec,i);
} // end of the sector digitization
for(i=0;iDelete();
delete row[i];
}
delete [] row; // delete the array of pointers to TObjArray-s
} // ntracks >0
} // end of Hits2DigitsSector
//_____________________________________________________________________________
void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
{
//-----------------------------------------------------------
// Single row digitization, coupling from the neighbouring
// rows taken into account
//-----------------------------------------------------------
//-----------------------------------------------------------------
// Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
// Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
//-----------------------------------------------------------------
Float_t zerosup = fTPCParam->GetZeroSup();
// Int_t nrows =fTPCParam->GetNRow(isec);
fCurrentIndex[1]= isec;
Int_t nofPads = fTPCParam->GetNPads(isec,irow);
Int_t nofTbins = fTPCParam->GetMaxTBin();
Int_t indexRange[4];
//
// Integrated signal for this row
// and a single track signal
//
AliTPCFastMatrix *m1 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // integrated
AliTPCFastMatrix *m2 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // single
//
AliTPCFastMatrix &total = *m1;
// Array of pointers to the label-signal list
Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
Float_t **pList = new Float_t* [nofDigits];
Int_t lp;
Int_t i1;
for(lp=0;lpGetNCrossRows(),0);
//Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
Int_t row1=irow;
Int_t row2=irow+2;
for (Int_t row= row1;row<=row2;row++){
Int_t nTracks= rows[row]->GetEntries();
for (i1=0;i1Zero(); // clear single track signal matrix
Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
GetList(trackLabel,nofPads,m2,indexRange,pList);
}
else GetSignal(rows[row],i1,0,m1,indexRange);
}
}
Int_t tracks[3];
AliDigits *dig = fDigitsArray->GetRow(isec,irow);
Int_t gi=-1;
Float_t fzerosup = zerosup+0.5;
for(Int_t it=0;it fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
}
else {
if(q <= 0.) continue; // do not fill zeros
if(q>2000.) q=2000.;
q *= 16.;
q = TMath::Nint(q);
}
//
// "real" signal or electronic noise (list = -1)?
//
for(Int_t j1=0;j1<3;j1++){
tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
}
//Begin_Html
/*
using of AliDigits object
*/
//End_Html
dig->SetDigitFast((Short_t)q,it,ip);
if (fDigitsArray->IsSimulated())
{
((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
}
} // end of loop over time buckets
} // end of lop over pads
//
// This row has been digitized, delete nonused stuff
//
for(lp=0;lpAt(ntr); // pointer to a track
AliTPCFastVector &v = *tv;
Float_t label = v(0);
Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1)-1)/2;
Int_t nElectrons = (tv->GetNrows()-1)/4;
indexRange[0]=9999; // min pad
indexRange[1]=-1; // max pad
indexRange[2]=9999; //min time
indexRange[3]=-1; // max time
AliTPCFastMatrix &signal = *m1;
AliTPCFastMatrix &total = *m2;
//
// Loop over all electrons
//
for(Int_t nel=0; nelGetTotalNormFac();
Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
Int_t *index = fTPCParam->GetResBin(0);
Float_t *weight = & (fTPCParam->GetResWeight(0));
if (n>0) for (Int_t i =0; i=0){
Int_t time=index[2];
Float_t qweight = *(weight)*eltoadcfac;
if (m1!=0) signal.UncheckedAt(pad,time)+=qweight;
total.UncheckedAt(pad,time)+=qweight;
if (indexRange[0]>pad) indexRange[0]=pad;
if (indexRange[1]time) indexRange[2]=time;
if (indexRange[3]