/************************************************************************** * 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. * **************************************************************************/ /* $Log$ Revision 1.8 2002/03/29 06:57:45 kowal2 Restored backward compatibility to use the hits from Dec. 2000 production. Revision 1.7 2001/10/21 19:04:55 hristov Several patches were done to adapt the barel reconstruction to the multi-event case. Some memory leaks were corrected. (Yu.Belikov) Revision 1.6 2001/08/30 09:28:48 hristov TTree names are explicitly set via SetName(name) and then Write() is called Revision 1.5 2001/07/20 14:32:44 kowal2 Processing of many events possible now Revision 1.4 2001/04/17 08:06:27 hristov Possibility to define the magnetic field in the reconstruction (Yu.Belikov) Revision 1.3 2000/10/05 16:14:01 kowal2 Forward declarations. Revision 1.2 2000/06/30 12:07:50 kowal2 Updated from the TPC-PreRelease branch Revision 1.1.2.1 2000/06/25 08:53:55 kowal2 Splitted from AliTPCtracking */ //------------------------------------------------------- // Implementation of the TPC clusterer // // Origin: Jouri Belikov, CERN, Jouri.Belikov@cern.ch //------------------------------------------------------- #include "AliTPCclusterer.h" #include "AliTPCcluster.h" #include #include #include "AliTPCClustersArray.h" #include "AliTPCClustersRow.h" #include "AliDigits.h" #include "AliSimDigits.h" #include "AliTPCParam.h" #include #include 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<IsOpen()) { cerr<<"AliTPC::Digits2Clusters(): output file not open !\n"; return; } const Int_t kMAXZ=par->GetMaxTBin()+2; char dname[100]; char cname[100]; if (eventn==-1) { // for backward compatibility sprintf(dname,"TreeD_75x40_100x60_150x60"); sprintf(cname,"TreeC_TPC"); } else { sprintf(dname,"TreeD_75x40_100x60_150x60_%d",eventn); sprintf(cname,"TreeC_TPC_%d",eventn); } TTree *t = (TTree *)gDirectory->Get(dname); if (!t) { cerr<<"Input tree with "<GetBranch("Segment")->SetAddress(&dummy); Stat_t nentries = t->GetEntries(); of->cd(); ((AliTPCParam*)par)->Write(par->GetTitle()); AliTPCClustersArray carray; carray.Setup(par); carray.SetClusterType("AliTPCcluster"); carray.MakeTree(); Int_t nclusters=0; for (Int_t n=0; nGetEvent(n); Int_t sec, row; if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) { cerr<<"AliTPC warning: invalid segment ID ! "<GetPadRowRadii(sec,row); Int_t npads, sign; { const Int_t kNIS=par->GetNInnerSector(), kNOS=par->GetNOuterSector(); if (sec < kNIS) { npads = par->GetNPadsLow(row); sign = (sec < kNIS/2) ? 1 : -1; } else { npads = par->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=par->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)*par->GetPadPitchWidth(sec)); c.SetZ(par->GetZWidth()*(c.GetZ()-1)); c.SetZ(c.GetZ() - 3.*par->GetZSigma()); // PASA delay c.SetZ(sign*(par->GetZLength() - 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++; } } carray.StoreRow(sec,row); carray.ClearRow(sec,row); //cerr<<"sector, row, compressed digits, clusters: " //<SetName(cname); carray.GetTree()->Write(); savedir->cd(); delete t; //Thanks to Mariana Bondila }