AliRoot Core  edcc906 (edcc906)
AliEMCALClusterizerv1.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
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  **************************************************************************/
15 
16 // --- ROOT system ---
17 #include <TFile.h>
18 #include <TMath.h>
19 #include <TMinuit.h>
20 #include <TTree.h>
21 #include <TBenchmark.h>
22 #include <TBrowser.h>
23 #include <TROOT.h>
24 #include <TList.h>
25 #include <TClonesArray.h>
26 
27 // --- Standard library ---
28 #include <cassert>
29 
30 // --- AliRoot header files ---
31 #include "AliLog.h"
32 #include "AliEMCALClusterizerv1.h"
33 #include "AliEMCALRecPoint.h"
34 #include "AliEMCALDigit.h"
35 #include "AliEMCALGeometry.h"
36 #include "AliCaloCalibPedestal.h"
37 #include "AliEMCALCalibData.h"
38 #include "AliEMCALCalibTime.h"
39 #include "AliESDCaloCluster.h"
40 #include "AliEMCALUnfolding.h"
41 
43 ClassImp(AliEMCALClusterizerv1) ;
45 
48 //____________________________________________________________________________
50 {
51  Init();
52 }
53 
58 //____________________________________________________________________________
60  : AliEMCALClusterizer(geometry)
61 { }
62 
70 //____________________________________________________________________________
72  AliEMCALCalibData * calib,
73  AliEMCALCalibTime * calibt,
74  AliCaloCalibPedestal * caloped)
75 : AliEMCALClusterizer(geometry, calib, calibt, caloped)
76 { }
77 
80 //____________________________________________________________________________
82 { }
83 
87 //____________________________________________________________________________
89 {
90  if(strstr(option,"tim"))
91  gBenchmark->Start("EMCALClusterizer");
92 
93  if(strstr(option,"print"))
94  Print("");
95 
96  //Get calibration parameters from file or digitizer default values.
98 
99  //Get dead channel map from file or digitizer default values.
101 
103 
104  MakeClusters(); //only the real clusters
105 
106  if(fToUnfold)
107  {
110  }
111 
112  //Evaluate position, dispersion and other RecPoint properties for EC section
113  Int_t index;
114  for(index = 0; index < fRecPoints->GetEntries(); index++)
115  {
116  AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index));
117 
118  if(rp)
119  {
121 
122  // Calculate the number of local maxima in cluster
123  // Do not do it for unfolded clusters
124  if(!fToUnfold)
125  {
127 
128  rp->SetNExMax(nMax);
129  }
130 
131  // For each rec.point set the distance to the nearest bad crystal
132  if (fCaloPed)
134  }
135  else AliFatal("Null rec point in list!");
136  }
137 
138  fRecPoints->Sort();
139 
140  for(index = 0; index < fRecPoints->GetEntries(); index++)
141  {
142  AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index));
143  if(rp)
144  {
145  rp->SetIndexInList(index);
146  rp->Print();
147  }
148  else AliFatal("Null rec point in list!");
149  }
150 
151  if (fTreeR)
152  fTreeR->Fill();
153 
154  if(strstr(option,"deb") || strstr(option,"all"))
155  PrintRecPoints(option);
156 
157  AliDebug(1,Form("EMCAL Clusterizer found %d Rec Points",fRecPoints->GetEntriesFast()));
158 
159  if(strstr(option,"tim"))
160  {
161  gBenchmark->Stop("EMCALClusterizer");
162  printf("Exec took %f seconds for Clusterizing",
163  gBenchmark->GetCpuTime("EMCALClusterizer"));
164  }
165 }
166 
175 //____________________________________________________________________________
176 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2, Bool_t & shared) const
177 {
178  Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
179  Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
180 
181  shared = kFALSE;
182 
183  fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
184  fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
185  fGeom->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, iphi1,ieta1);
186  fGeom->GetCellPhiEtaIndexInSModule(nSupMod2,nModule2,nIphi2,nIeta2, iphi2,ieta2);
187 
188  //If different SM, check if they are in the same phi, then consider cells close to eta=0 as neighbours; May 2010
189  if (nSupMod1 != nSupMod2 )
190  {
191  //Check if the 2 SM are in the same PHI position (0,1), (2,3), ...
192  Float_t smPhi1 = fGeom->GetEMCGeometry()->GetPhiCenterOfSM(nSupMod1);
193  Float_t smPhi2 = fGeom->GetEMCGeometry()->GetPhiCenterOfSM(nSupMod2);
194 
195  if(!TMath::AreEqualAbs(smPhi1, smPhi2, 1e-3)) return 2; //Not phi rack equal, not neighbours
196 
197  // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
198  // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
199  if(nSupMod1%2) ieta1+=AliEMCALGeoParams::fgkEMCALCols;
201 
202  shared = kTRUE; // maybe a shared cluster, we know this later, set it for the moment.
203  } //Different SM, same phi
204 
205  Int_t rowdiff = TMath::Abs(iphi1 - iphi2);
206  Int_t coldiff = TMath::Abs(ieta1 - ieta2);
207 
208  // neighbours with at least common side; May 11, 2007
209  if ((coldiff==0 && TMath::Abs(rowdiff)==1) || (rowdiff==0 && TMath::Abs(coldiff)==1))
210  {
211  //Diagonal?
212  //if ((coldiff==0 && TMath::Abs(rowdiff==1)) || (rowdiff==0 && TMath::Abs(coldiff==1)) || (TMath::Abs(rowdiff)==1 && TMath::Abs(coldiff==1))) rv = 1;
213 
214  if (gDebug == 2)
215  printf("AliEMCALClusterizerv1::AreNeighbours(): id1=%d, (row %d, col %d) ; id2=%d, (row %d, col %d), shared %d \n",
216  d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2, shared);
217  return 1;
218  } //Neighbours
219  else
220  {
221  shared = kFALSE;
222  return 2;
223  } //Not neighbours
224 }
225 
229 //____________________________________________________________________________
231 {
232  if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
233 
234  fRecPoints->Delete();
235 
236  // Set up TObjArray with pointers to digits to work on calibrated digits
237  TObjArray *digitsC = new TObjArray();
238  AliEMCALDigit *digit;
239  Float_t dEnergyCalibrated = 0.0, ehs = 0.0, time = 0.0;
240 
241  TIter nextdigit(fDigitsArr);
242  while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) )
243  {
244  // Calibrate and clean up digits
245  dEnergyCalibrated = digit->GetAmplitude();
246  time = digit->GetTime();
247  Calibrate(dEnergyCalibrated,time,digit->GetId());
248  digit->SetCalibAmp(dEnergyCalibrated);
249  digit->SetTime(time);
250 
251  if ( dEnergyCalibrated < fMinECut || time > fTimeMax || time < fTimeMin )
252  {
253  continue;
254  }
255  else if (!fGeom->CheckAbsCellId(digit->GetId()))
256  continue;
257  else
258  {
259  ehs += dEnergyCalibrated;
260  digitsC->AddLast(digit);
261  }
262  }
263 
264  AliDebug(1,Form("MakeClusters: Number of digits %d -> (e %f), ehs %f\n",
265  fDigitsArr->GetEntries(),fMinECut,ehs));
266 
267  TIter nextdigitC(digitsC);
268  while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) )
269  {
270  // Scan over the list of digitsC
271  TArrayI clusterECAdigitslist(fDigitsArr->GetEntries());
272  dEnergyCalibrated = digit->GetCalibAmp();
273  time = digit->GetTime();
274 
275  if(fGeom->CheckAbsCellId(digit->GetId()) && ( dEnergyCalibrated > fECAClusteringThreshold ) )
276  {
277  // Start a new Tower RecPoint
278  if(fNumberOfECAClusters >= fRecPoints->GetSize()) fRecPoints->Expand(2*fNumberOfECAClusters+1);
279 
280  AliEMCALRecPoint *recPoint = new AliEMCALRecPoint("");
281  fRecPoints->AddAt(recPoint, fNumberOfECAClusters);
282  recPoint = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(fNumberOfECAClusters));
283 
284  if (recPoint)
285  {
287 
289  recPoint->AddDigit(*digit, digit->GetCalibAmp(), kFALSE); //Time or TimeR?
290  TObjArray clusterDigits;
291  clusterDigits.AddLast(digit);
292  digitsC->Remove(digit);
293 
294  AliDebug(1,Form("MakeClusters: OK id = %d, ene = %f , cell.th. = %f \n", digit->GetId(), dEnergyCalibrated, fECAClusteringThreshold)); //Time or TimeR?
295 
296  // Grow cluster by finding neighbours
297  TIter nextClusterDigit(&clusterDigits);
298 
299  while ( (digit = dynamic_cast<AliEMCALDigit*>(nextClusterDigit())) )
300  {
301  // Scan over digits in cluster
302  TIter nextdigitN(digitsC);
303  AliEMCALDigit *digitN = 0; // digi neighbor
304  while ( (digitN = (AliEMCALDigit *)nextdigitN()) )
305  {
306  // Scan over all digits to look for neighbours
307  // Do not add digits with too different time
308  Bool_t shared = kFALSE;//cluster shared by 2 SuperModules?
309  if(TMath::Abs(time - digitN->GetTime()) > fTimeCut ) continue; //Time or TimeR?
310 
311  if (AreNeighbours(digit, digitN, shared)==1)
312  { // call (digit,digitN) in THAT order !!!!!
313  recPoint->AddDigit(*digitN, digitN->GetCalibAmp(), shared);
314  clusterDigits.AddLast(digitN);
315  digitsC->Remove(digitN);
316  } // if(ineb==1)
317  } // scan over digits
318  } // scan over digits already in cluster
319 
320  AliDebug(2,Form("MakeClusters: %d digitd, energy %f \n", clusterDigits.GetEntries(), recPoint->GetEnergy()));
321  }//recpoint
322  else AliFatal("Null recpoint in array!");
323  } // If seed found
324  } // while digit
325 
326  delete digitsC;
327 
328  AliDebug(1,Form("total no of clusters %d from %d digits",fNumberOfECAClusters,fDigitsArr->GetEntriesFast()));
329 }
Bool_t CheckAbsCellId(Int_t absId) const
void EvalDistanceToBadChannels(AliCaloCalibPedestal *caloped)
printf("Chi2/npoints = %f\n", TMath::Sqrt(chi2/npoints))
virtual Int_t GetNumberOfLocalMax(Int_t nDigitMult, Float_t locMaxCut, TClonesArray *digits) const
Int_t fNumberOfECAClusters
number of clusters found in EC section
virtual void Calibrate(Float_t &amp, Float_t &time, const Int_t cellId)
virtual void Digits2Clusters(Option_t *option)
#define TObjArray
Float_t GetTime(void) const
Definition: AliEMCALDigit.h:64
TObjArray * fRecPoints
array with EMCAL clusters
void SetCalibAmp(Float_t amp)
Float_t fECAW0
logarithmic weight for the cluster center of gravity calculation
TBenchmark * gBenchmark
Float_t fECALocMaxCut
minimum energy difference to distinguish local maxima in a cluster
virtual void SetInput(Int_t numberOfECAClusters, TObjArray *recPoints, TClonesArray *digitsArr)
virtual Float_t GetEnergy() const
Cell energy calibration factors container class.
EMCal digits object.
Definition: AliEMCALDigit.h:30
Double_t GetPhiCenterOfSM(Int_t nsupmod) const
AliEMCALClusterizerv1()
Default constructor.
void GetCellPhiEtaIndexInSModule(Int_t nSupMod, Int_t nModule, Int_t nIphi, Int_t nIeta, Int_t &iphi, Int_t &ieta) const
void SetTime(Float_t time)
Definition: AliEMCALDigit.h:74
Float_t fTimeCut
maximum time difference between the digits inside EMC cluster
Int_t GetMultiplicity(void) const
Bool_t fJustClusters
false for standard reco
virtual void GetCalibrationParameters(void)
Bool_t fToUnfold
says if unfolding should be performed
Bool_t GetCellIndex(Int_t absId, Int_t &nSupMod, Int_t &nModule, Int_t &nIphi, Int_t &nIeta) const
virtual void MakeUnfolding()
virtual void AddDigit(AliEMCALDigit &digit, const Float_t energy, const Bool_t shared)
Float_t fTimeMin
minimum time of physical signal in a cell/digit
EMCal rec_points object.
AliEMCALGeometry * fGeom
! pointer to geometry for utilities
Cell time shifts container class.
Double_t GetCalibAmp() const
Clusterize neighbour cells, no split, unfolding possible.
pedestal/bad map monitoring and calibration tools
Float_t fTimeMax
maximum time of physical signal in a cell/digit
#define AliFatal(message)
Definition: AliLog.h:640
static const int fgkEMCALCols
Number of columns per module for EMCAL.
TClonesArray * fDigitsArr
array with EMCAL digits
#define AliDebug(logLevel, message)
Definition: AliLog.h:300
AliCaloCalibPedestal * fCaloPed
! tower status map if aval
virtual void Print(Option_t *option="") const
Print the list of digits belonging to the cluster.
AliEMCALEMCGeometry * GetEMCGeometry() const
Int_t GetId() const
Definition: AliDigitNew.h:23
virtual void GetCaloCalibPedestal(void)
AliEMCALUnfolding * fClusterUnfolding
! pointer to unfolding object
virtual void EvalAll(Float_t logWeight, TClonesArray *digits, const Bool_t justClusters)
Evaluates cluster parameters: position, shower shape, primaries ...
Float_t fECAClusteringThreshold
minimum energy to seed a EC digit in a cluster
TTree * fTreeR
tree with output clusters
virtual ~AliEMCALClusterizerv1()
Destructtor.
virtual void SetClusterType(Int_t ver)
void SetNExMax(Int_t nmax=1)
void SetIndexInList(Int_t val)
Float_t GetAmplitude() const
Definition: AliEMCALDigit.h:55
virtual void Print(Option_t *option) const
Print clusterizer parameters.
EMCal geometry, singleton.
virtual Int_t AreNeighbours(AliEMCALDigit *d1, AliEMCALDigit *d2, Bool_t &shared) const
Float_t fMinECut
minimum energy for a digit to be a member of a cluster
Base class for the clusterization algorithm (pure abstract)
virtual void PrintRecPoints(Option_t *option)