/************************************************************************** * 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$ */ //------------------------------------------------------- // Implementation of the TPC clusterer // // Origin: Jouri Belikov, CERN, Jouri.Belikov@cern.ch //------------------------------------------------------- #include "AliTPCclusterer.h" #include "AliTPCcluster.h" #include #include #include "AliDigits.h" #include "AliSimDigits.h" #include "AliTPCParam.h" #include "AliTPCClustersRow.h" #include ClassImp(AliTPCclusterer) //____________________________________________________ AliTPCclusterer::AliTPCclusterer(const AliTPCclusterer ¶m) :TObject(), fPar(0) { // // dummy // fPar = param.fPar; } //----------------------------------------------------- AliTPCclusterer & AliTPCclusterer::operator =(const AliTPCclusterer & param) { // // assignment operator - dummy // fPar = param.fPar; return (*this); } //____________________________________________________ void AliTPCclusterer::FindPeaks(Int_t k,Int_t max, AliBin *b,Int_t *idx,UInt_t *msk,Int_t& n) { //find local maxima if (n<31) if (IsMaximum(k,max,b)) { idx[n]=k; msk[n]=(2<GetBranch("Segment"); if (!branch) { Error("Digits2Cluster","Can't get the branch !"); return 1; } AliSimDigits digarr, *dummy=&digarr; branch->SetAddress(&dummy); AliTPCClustersRow ddd,*clrow=&ddd; clrow = new AliTPCClustersRow("AliTPCcluster"); cTree->Branch("Segment","AliTPCClustersRow",&clrow,32000,200); const Int_t kMAXZ=fPar->GetMaxTBin()+2; Int_t nclusters=0; Int_t nentries = (Int_t)dTree->GetEntries(); for (Int_t n=0; nGetEvent(n); if (!fPar->AdjustSectorRow(digarr.GetID(),sec,row)) { Error("Digits2Clusters","!nvalid segment ID ! %d",digarr.GetID()); continue; } clrow->SetID(digarr.GetID()); cTree->GetBranch("Segment")->SetAddress(&clrow); Float_t rx=fPar->GetPadRowRadii(sec,row); Int_t npads, sign; { const Int_t kNIS=fPar->GetNInnerSector(), kNOS=fPar->GetNOuterSector(); if (sec < kNIS) { npads = fPar->GetNPadsLow(row); sign = (sec < kNIS/2) ? 1 : -1; } else { npads = fPar->GetNPadsUp(row); sign = ((sec-kNIS) < kNOS/2) ? 1 : -1; } } const Int_t kMAXBIN=kMAXZ*(npads+2); AliBin *bins=new AliBin[kMAXBIN]; for (Int_t ii=0;iiGetZeroSup()) continue; Int_t j=digarr.CurrentRow()+1, i=digarr.CurrentColumn()+1; bins[i*kMAXZ+j].SetQ(dig); bins[i*kMAXZ+j].SetMask(1); } while (digarr.Next()); Int_t ncl=0; for (Int_t i=0; i30) continue; Int_t k,l; for (k=0; k1 || dj>1) continue; if (bins[idx[k]].GetQ() > bins[idx[l]].GetQ()) { msk[l]=msk[k]; idx[l]*=-1; } else { msk[k]=msk[l]; idx[k]*=-1; break; } } } for (k=0; kGetPadPitchWidth(sec); c.SetSigmaY2((s2 + 1./12.)*w*w); if (s2 != 0.) { c.SetSigmaY2(c.GetSigmaY2()*0.108); if (secGetNInnerSector()) c.SetSigmaY2(c.GetSigmaY2()*2.07); } s2 = c.GetSigmaZ2()/c.GetQ() - c.GetZ()*c.GetZ(); w=fPar->GetZWidth(); c.SetSigmaZ2((s2 + 1./12.)*w*w); if (s2 != 0.) { c.SetSigmaZ2(c.GetSigmaZ2()*0.169); if (secGetNInnerSector()) c.SetSigmaZ2(c.GetSigmaZ2()*1.77); } c.SetY((c.GetY() - 0.5 - 0.5*npads)*fPar->GetPadPitchWidth(sec)); c.SetZ(fPar->GetZWidth()*(c.GetZ()-1)); c.SetZ(c.GetZ() - 3.*fPar->GetZSigma()); // PASA delay c.SetZ(sign*(fPar->GetZLength(sec) - c.GetZ())); if (rx<230./250.*TMath::Abs(c.GetZ())) continue; Int_t ki=idx[k]/kMAXZ, kj=idx[k] - ki*kMAXZ; c.SetLabel(digarr.GetTrackID(kj-1,ki-1,0),0); c.SetLabel(digarr.GetTrackID(kj-1,ki-1,1),1); c.SetLabel(digarr.GetTrackID(kj-1,ki-1,2),2); c.SetQ(bins[idx[k]].GetQ()); if (ki==1 || ki==npads || kj==1 || kj==kMAXZ-2) { c.SetSigmaY2(c.GetSigmaY2()*25.); c.SetSigmaZ2(c.GetSigmaZ2()*4.); } clrow->InsertCluster(&c); ncl++; } } cTree->Fill(); clrow->GetArray()->Clear(); nclusters+=ncl; delete[] bins; } Info("Digits2Cluster","Number of found clusters : %d",nclusters); return 0; }