CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes | Friends
HFClusterAlgo Class Reference

#include <HFClusterAlgo.h>

Classes

struct  HFCompleteHit
 

Public Member Functions

void clusterize (const HFRecHitCollection &hf, const CaloGeometry &geom, reco::HFEMClusterShapeCollection &clusters, reco::SuperClusterCollection &SuperClusters)
 
 HFClusterAlgo ()
 
void isMC (bool isMC)
 
void resetForRun ()
 
void setup (double minTowerEnergy, double seedThreshold, double maximumSL, double m_maximumRenergy, bool usePMTflag, bool usePulseflag, bool forcePulseFlagMC, int correctionSet)
 

Private Member Functions

bool isPMTHit (const HFRecHit &hfr)
 
bool makeCluster (const HcalDetId &seedid, const HFRecHitCollection &hf, const CaloGeometry &geom, reco::HFEMClusterShape &clusShp, reco::SuperCluster &SClus)
 

Private Attributes

std::vector< double > m_correctionByEta
 
int m_correctionSet
 
std::vector< double > m_cutByEta
 
bool m_forcePulseFlagMC
 
bool m_isMC
 
double m_maximumRenergy
 
double m_maximumSL
 
double m_minTowerEnergy
 
std::vector< double > m_seedmnEta
 
std::vector< double > m_seedmnPhi
 
std::vector< double > m_seedMXeta
 
std::vector< double > m_seedMXphi
 
double m_seedThreshold
 
bool m_usePMTFlag
 
bool m_usePulseFlag
 

Friends

class CompareHFCompleteHitET
 
class CompareHFCore
 

Detailed Description

Author
Kevin Klapoetke (Minnesota)

$Id:version 1.2

Author
K. Klapoetke – Minnesota

Definition at line 22 of file HFClusterAlgo.h.

Constructor & Destructor Documentation

HFClusterAlgo::HFClusterAlgo ( )

Definition at line 20 of file HFClusterAlgo.cc.

21 {
22  m_isMC=true; // safest
23 
24 }

Member Function Documentation

void HFClusterAlgo::clusterize ( const HFRecHitCollection hf,
const CaloGeometry geom,
reco::HFEMClusterShapeCollection clusters,
reco::SuperClusterCollection SuperClusters 
)

Analyze the hits

Definition at line 91 of file HFClusterAlgo.cc.

References edm::SortedCollection< T, SORT >::begin(), edm::SortedCollection< T, SORT >::end(), HFClusterAlgo::HFCompleteHit::energy, HFClusterAlgo::HFCompleteHit::et, eta(), edm::SortedCollection< T, SORT >::find(), relativeConstraints::geom, CaloCellGeometry::getCorners(), CaloGeometry::getGeometry(), CaloGeometry::getPosition(), HcalForward, i, HFClusterAlgo::HFCompleteHit::id, indexByEta(), j, relval_steps::k, EZArrayFL< T >::size(), python.multivaluedict::sort(), and detailsBasic3DVector::z.

Referenced by HFEMClusterProducer::produce().

94  {
95 
96  std::vector<HFCompleteHit> protoseeds, seeds;
98  std::vector<HFCompleteHit>::iterator i;
99  std::vector<HFCompleteHit>::iterator k;
100  int dP, dE, PWrap;
101  bool isok=true;
102  HFEMClusterShape clusShp;
103 
104  SuperCluster Sclus;
105  bool doCluster=false;
106 
107  for (j=hf.begin(); j!= hf.end(); j++) {
108  const int aieta=j->id().ietaAbs();
109  int iz=(aieta-29);
110  // only long fibers and not 29,40,41 allowed to be considered as seeds
111  if (j->id().depth()!=1) continue;
112  if (aieta==40 || aieta==41 || aieta==29) continue;
113 
114 
115  if (iz<0 || iz>12) {
116  edm::LogWarning("HFClusterAlgo") << "Strange invalid HF hit: " << j->id();
117  continue;
118  }
119 
120  if (m_cutByEta[iz]<0) {
121  double eta=geom.getPosition(j->id()).eta();
122  m_cutByEta[iz]=m_seedThreshold*cosh(eta); // convert ET to E for this ring
123  const CaloCellGeometry* ccg=geom.getGeometry(j->id());
124  const CaloCellGeometry::CornersVec& CellCorners=ccg->getCorners();
125  for(size_t sc=0;sc<CellCorners.size();sc++){
126  if(fabs(CellCorners[sc].z())<1200){
127  if(fabs(CellCorners[sc].eta())<m_seedmnEta[iz])m_seedmnEta[iz]=fabs(CellCorners[sc].eta());
128  if(fabs(CellCorners[sc].eta())>m_seedMXeta[iz])m_seedMXeta[iz]=fabs(CellCorners[sc].eta());
129  }
130  }
131  }
132  double elong=j->energy()*m_correctionByEta[indexByEta(j->id())];
133  if (elong>m_cutByEta[iz]) {
134  j2=hf.find(HcalDetId(HcalForward,j->id().ieta(),j->id().iphi(),2));
135  double eshort=(j2==hf.end())?(0):(j2->energy());
136  if (j2!=hf.end())
137  eshort*=m_correctionByEta[indexByEta(j2->id())];
138  if (((elong-eshort)/(elong+eshort))>m_maximumSL) continue;
139  //if ((m_usePMTFlag)&&(j->flagField(4,1))) continue;
140  //if ((m_usePulseFlag)&&(j->flagField(1,1))) continue;
141  if(isPMTHit(*j)) continue;
142 
143  HFCompleteHit ahit;
144  double eta=geom.getPosition(j->id()).eta();
145  ahit.id=j->id();
146  ahit.energy=elong;
147  ahit.et=ahit.energy/cosh(eta);
148  protoseeds.push_back(ahit);
149  }
150  }
151 
152  if(!protoseeds.empty()){
153  std::sort(protoseeds.begin(), protoseeds.end(), CompareHFCompleteHitET());
154  for (i=protoseeds.begin(); i!= protoseeds.end(); i++) {
155  isok=true;
156  doCluster=false;
157 
158  if ( (i==protoseeds.begin()) && (isok) ) {
159  doCluster=true;
160  }else {
161  // check for overlap with existing clusters
162  for (k=seeds.begin(); isok && k!=seeds.end(); k++) { //i->hits, k->seeds
163 
164  for (dE=-2; dE<=2; dE++)
165  for (dP=-4;dP<=4; dP+=2) {
166  PWrap=k->id.iphi()+dP;
167  if (PWrap<0)
168  PWrap+=72;
169  if (PWrap>72)
170  PWrap-=72;
171 
172  if ( (i->id.iphi()==PWrap) && (i->id.ieta()==k->id.ieta()+dE))
173  isok = false;
174  }
175  }
176  if (isok) {
177  doCluster=true;
178  }
179  }
180  if (doCluster) {
181  seeds.push_back(*i);
182 
183  bool clusterOk=makeCluster( i->id(),hf, geom,clusShp,Sclus);
184  if (clusterOk) { // cluster is _not_ ok if seed is rejected due to other cuts
185  clusterShapes.push_back(clusShp);
186  SuperClusters.push_back(Sclus);
187  }
188 
189  }
190  }//end protoseed loop
191  }//end if seeCount
192 }
int i
Definition: DBlmapReader.cc:9
std::vector< double > m_cutByEta
Definition: HFClusterAlgo.h:48
double m_maximumSL
Definition: HFClusterAlgo.h:43
std::vector< double > m_correctionByEta
Definition: HFClusterAlgo.h:49
std::vector< HFRecHit >::const_iterator const_iterator
std::vector< double > m_seedmnEta
Definition: HFClusterAlgo.h:50
friend class CompareHFCompleteHitET
Definition: HFClusterAlgo.h:40
T eta() const
float float float z
std::vector< double > m_seedMXeta
Definition: HFClusterAlgo.h:51
bool makeCluster(const HcalDetId &seedid, const HFRecHitCollection &hf, const CaloGeometry &geom, reco::HFEMClusterShape &clusShp, reco::SuperCluster &SClus)
int j
Definition: DBlmapReader.cc:9
const GlobalPoint & getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:68
const_iterator end() const
static int indexByEta(HcalDetId id)
size_type size() const
Definition: EZArrayFL.h:81
iterator find(key_type k)
const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:76
double m_seedThreshold
Definition: HFClusterAlgo.h:43
const CornersVec & getCorners() const
Returns the corner points of this cell&#39;s volume.
bool isPMTHit(const HFRecHit &hfr)
const_iterator begin() const
void HFClusterAlgo::isMC ( bool  isMC)
inline

Definition at line 28 of file HFClusterAlgo.h.

References isMC(), and m_isMC.

Referenced by isMC(), and HFEMClusterProducer::produce().

28 { m_isMC=isMC; }
void isMC(bool isMC)
Definition: HFClusterAlgo.h:28
bool HFClusterAlgo::isPMTHit ( const HFRecHit hfr)
private

Definition at line 398 of file HFClusterAlgo.cc.

References CaloRecHit::flagField(), HcalCaloFlagLabels::HFDigiTime, and HcalCaloFlagLabels::HFLongShort.

398  {
399 
400  bool pmthit=false;
401 
402  if((hfr.flagField(HcalCaloFlagLabels::HFLongShort))&&(m_usePMTFlag)) pmthit=true;
403  if (!(m_isMC && !m_forcePulseFlagMC))
404  if((hfr.flagField(HcalCaloFlagLabels::HFDigiTime))&&(m_usePulseFlag)) pmthit=true;
405 
406  return pmthit;
407 }
bool m_forcePulseFlagMC
Definition: HFClusterAlgo.h:45
uint32_t flagField(int base, int width=1) const
Definition: CaloRecHit.cc:31
bool HFClusterAlgo::makeCluster ( const HcalDetId seedid,
const HFRecHitCollection hf,
const CaloGeometry geom,
reco::HFEMClusterShape clusShp,
reco::SuperCluster SClus 
)
private

Definition at line 195 of file HFClusterAlgo.cc.

References funct::abs(), edm::SortedCollection< T, SORT >::end(), PV3DBase< T, PVType, FrameType >::eta(), eta(), edm::SortedCollection< T, SORT >::find(), CaloGeometry::getPosition(), HcalForward, i, customizeTrackingMonitorSeedNumber::idx, HcalDetId::ieta(), HcalDetId::ietaAbs(), indexByEta(), HcalDetId::iphi(), create_public_lumi_plots::log, M_PI, bookConverter::max, AlCaHLTBitMon_ParallelJobs::p, phi, PV3DBase< T, PVType, FrameType >::phi(), python.multivaluedict::sort(), w, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and HcalDetId::zside().

199  {
200 
201 
202  double w=0;//sum over all log E's
203  double wgt=0;
204  double w_e=0;//sum over ieat*energy
205  double w_x=0;
206  double w_y=0;
207  double w_z=0;
208  double wp_e=0;//sum over iphi*energy
209  double e_e=0;//nonwieghted eta sum
210  double e_ep=0; //nonweighted phi sum
211 
212  double l_3=0;//sum for enenergy in 3x3 long fibers etc.
213  double s_3=0;
214  double l_5=0;
215  double s_5=0;
216  double l_1=0;
217  double s_1=0;
218  int de, dp, phiWrap;
219  double l_1e=0;
220  GlobalPoint sp=geom.getPosition(seedid);
221  std::vector<double> coreCanid;
222  std::vector<double>::const_iterator ci;
224  std::vector<DetId> usedHits;
225 
227  HcalDetId sid(HcalForward,seedid.ieta(),seedid.iphi(),1);
228  si=hf.find(sid);
229 
230  bool clusterOk=true; // assume the best to start...
231 
232  // lots happens here
233  // edge type 1 has 40/41 in 3x3 and 5x5
234  bool edge_type1=seedid.ietaAbs()==39 && (seedid.iphi()%4)==3;
235 
236  double e_seed=si->energy()*m_correctionByEta[indexByEta(si->id())];
237 
238  for (de=-2; de<=2; de++)
239  for (dp=-4;dp<=4; dp+=2) {
240  phiWrap=seedid.iphi()+dp;
241  if (phiWrap<0)
242  phiWrap+=72;
243  if (phiWrap>72)
244  phiWrap-=72;
245 
246 
247  /* Handling of phi-width change problems */
248  if (edge_type1 && de==seedid.zside()) {
249  if (dp==-2) { // we want it in the 3x3
250  phiWrap-=2;
251  if (phiWrap<0)
252  phiWrap+=72;
253  }
254  else if (dp==-4) {
255  continue; // but not double counted in 5x5
256  }
257  }
258 
259  HcalDetId idl(HcalForward,seedid.ieta()+de,phiWrap,1);
260  HcalDetId ids(HcalForward,seedid.ieta()+de,phiWrap,2);
261 
262 
263  il=hf.find(idl);
264  is=hf.find(ids);
265 
266 
267 
268 
269  double e_long=1.0;
270  double e_short=0.0;
271  if (il!=hf.end()) e_long=il->energy()*m_correctionByEta[indexByEta(il->id())];
272  if (e_long <= m_minTowerEnergy) e_long=0.0;
273  if (is!=hf.end()) e_short=is->energy()*m_correctionByEta[indexByEta(is->id())];
274  if (e_short <= m_minTowerEnergy) e_short=0.0;
275  double eRatio=(e_long-e_short)/std::max(1.0,(e_long+e_short));
276 
277  // require S/L > a minimum amount for inclusion
278  if ((abs(eRatio) > m_maximumSL)&&(std::max(e_long,e_short) > m_maximumRenergy)) {
279  if (dp==0 && de==0) clusterOk=false; // somehow, the seed is hosed
280  continue;
281  }
282 
283  if((il!=hf.end())&&(isPMTHit(*il))){
284  if (dp==0 && de==0) clusterOk=false; // somehow, the seed is hosed
285  continue;//continue to next hit, do not include this one in cluster
286  }
287 
288 
289 
290 
291 
292  if (e_long > m_minTowerEnergy && il!=hf.end()) {
293 
294  // record usage
295  usedHits.push_back(idl.rawId());
296  // always in the 5x5
297  l_5+=e_long;
298  // maybe in the 3x3
299  if ((de>-2)&&(de<2)&&(dp>-4)&&(dp<4)) {
300  l_3+=e_long;
301 
302  // sometimes in the 1x1
303  if ((dp==0)&&(de==0)) {
304  l_1=e_long;
305  }
306 
307  // maybe in the core?
308  if ((de>-2)&&(de<2)&&(dp>-4)&&(dp<4)&&(e_long>(.5*e_seed))) {
309  coreCanid.push_back(e_long);
310  }
311 
312  // position calculation
313  GlobalPoint p=geom.getPosition(idl);
314 
315  double d_p = p.phi()-sp.phi();
316  while (d_p < -M_PI)
317  d_p+=2*M_PI;
318  while (d_p > M_PI)
319  d_p-=2*M_PI;
320  double d_e = p.eta()-sp.eta();
321 
322  wgt=log((e_long));
323  if (wgt>0){
324  w+=wgt;
325  w_e+=(d_e)*wgt;
326  wp_e+=(d_p)*wgt;
327  e_e+=d_e;
328  e_ep+=d_p;
329 
330  w_x+=(p.x())*wgt;//(p.x()-sp.x())*wgt;
331  w_y+=(p.y())*wgt;
332  w_z+=(p.z())*wgt;
333  }
334  }
335  } else {
336  if (dp==0 && de==0) clusterOk=false; // somehow, the seed is hosed
337  }
338 
339  if (e_short > m_minTowerEnergy && is!=hf.end()) {
340  // record usage
341  usedHits.push_back(ids.rawId());
342  // always in the 5x5
343  s_5+=e_short;
344  // maybe in the 3x3
345  if ((de>-2)&&(de<2)&&(dp>-4)&&(dp<4)) {
346  s_3+=e_short;
347  }
348  // sometimes in the 1x1
349  if ((dp==0)&&(de==0)) {
350  s_1=e_short;
351  }
352  }
353  }
354 
355 
356  if (!clusterOk) return false;
357 
358  //Core sorting done here
359  std::sort(coreCanid.begin(), coreCanid.end(), CompareHFCore());
360  for (ci=coreCanid.begin();ci!=coreCanid.end();ci++){
361  if(ci==coreCanid.begin()){
362  l_1e=*ci;
363  }else if (*ci>.5*l_1e){
364  l_1e+=*ci;
365  }
366  }//core sorting end
367 
368  double z_=w_z/w;
369  double x_=w_x/w;
370  double y_=w_y/w;
371  math::XYZPoint xyzclus(x_,y_,z_);
372  //calcualte position, final
373  double eta=xyzclus.eta();//w_e/w+sp.eta();
374 
375  double phi=xyzclus.phi();//(wp_e/w)+sp.phi();
376 
377  while (phi < -M_PI)
378  phi+=2*M_PI;
379  while (phi > M_PI)
380  phi-=2*M_PI;
381 
382  //calculate cell phi and cell eta
383  int idx= fabs(seedid.ieta())-29;
384  int ipx=seedid.iphi();
385  double Cphi =(phi-m_seedmnPhi[ipx])/(m_seedMXphi[ipx]-m_seedmnPhi[ipx]);
386  double Ceta=(fabs(eta)- m_seedmnEta[idx])/(m_seedMXeta[idx]-m_seedmnEta[idx]);
387 
388  //return HFEMClusterShape, SuperCluster
389  HFEMClusterShape myClusShp(l_1, s_1, l_3, s_3, l_5,s_5, l_1e,Ceta, Cphi,seedid);
390  clusShp = myClusShp;
391 
392  SuperCluster MySclus(l_3,xyzclus);
393  Sclus=MySclus;
394 
395  return clusterOk;
396 
397 }
double m_minTowerEnergy
Definition: HFClusterAlgo.h:43
int i
Definition: DBlmapReader.cc:9
double m_maximumRenergy
Definition: HFClusterAlgo.h:43
const double w
Definition: UKUtility.cc:23
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.h:32
double m_maximumSL
Definition: HFClusterAlgo.h:43
std::vector< double > m_correctionByEta
Definition: HFClusterAlgo.h:49
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
std::vector< HFRecHit >::const_iterator const_iterator
T y() const
Definition: PV3DBase.h:63
std::vector< double > m_seedmnEta
Definition: HFClusterAlgo.h:50
T eta() const
friend class CompareHFCore
Definition: HFClusterAlgo.h:41
std::vector< double > m_seedMXeta
Definition: HFClusterAlgo.h:51
T z() const
Definition: PV3DBase.h:64
int ieta() const
get the cell ieta
Definition: HcalDetId.h:36
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const GlobalPoint & getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:68
std::vector< double > m_seedmnPhi
Definition: HFClusterAlgo.h:52
#define M_PI
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:34
const_iterator end() const
static int indexByEta(HcalDetId id)
int iphi() const
get the cell iphi
Definition: HcalDetId.h:38
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
std::vector< double > m_seedMXphi
Definition: HFClusterAlgo.h:53
T eta() const
Definition: PV3DBase.h:76
iterator find(key_type k)
T x() const
Definition: PV3DBase.h:62
bool isPMTHit(const HFRecHit &hfr)
Definition: DDAxes.h:10
void HFClusterAlgo::resetForRun ( )

Definition at line 408 of file HFClusterAlgo.cc.

References cuy::ii.

Referenced by HFEMClusterProducer::beginRun().

408  {
409  edm::LogInfo("HFClusterAlgo")<<"Resetting for Run!";
410  for(int ii=0;ii<13;ii++){
411  m_cutByEta.push_back(-1);
412  m_seedmnEta.push_back(99);
413  m_seedMXeta.push_back(-1);
414  }
415 }
std::vector< double > m_cutByEta
Definition: HFClusterAlgo.h:48
std::vector< double > m_seedmnEta
Definition: HFClusterAlgo.h:50
int ii
Definition: cuy.py:588
std::vector< double > m_seedMXeta
Definition: HFClusterAlgo.h:51
void HFClusterAlgo::setup ( double  minTowerEnergy,
double  seedThreshold,
double  maximumSL,
double  m_maximumRenergy,
bool  usePMTflag,
bool  usePulseflag,
bool  forcePulseFlagMC,
int  correctionSet 
)

Definition at line 48 of file HFClusterAlgo.cc.

References cuy::ii, M_PI, and MCMaterialCorrections_3XX.

Referenced by HFEMClusterProducer::HFEMClusterProducer().

49  {
50  m_seedThreshold=seedThreshold;
51  m_minTowerEnergy=minTowerEnergy;
52  m_maximumSL=maximumSL;
53  m_usePMTFlag=usePMTflag;
54  m_usePulseFlag=usePulseflag;
55  m_forcePulseFlagMC=forcePulseFlagMC;
56  m_maximumRenergy=maximumRenergy;
57  m_correctionSet = correctionSet;
58 
59  for(int ii=0;ii<13;ii++){
60  m_cutByEta.push_back(-1);
61  m_seedmnEta.push_back(99);
62  m_seedMXeta.push_back(-1);
63  }
64  for(int ii=0;ii<73;ii++){
65  double minphi=0.0872664*(ii-1);
66  double maxphi=0.0872664*(ii+1);
67  while (minphi < -M_PI)
68  minphi+=2*M_PI;
69  while (minphi > M_PI)
70  minphi-=2*M_PI;
71  while (maxphi < -M_PI)
72  maxphi+=2*M_PI;
73  while (maxphi > M_PI)
74  maxphi-=2*M_PI;
75  if(ii==37)
76  minphi=-3.1415904;
77  m_seedmnPhi.push_back(minphi);
78  m_seedMXphi.push_back(maxphi);
79  }
80 
81  // always set all the corrections to one...
82  for (int ii=0; ii<13*2; ii++)
83  m_correctionByEta.push_back(1.0);
84  if (m_correctionSet==1) { // corrections for material from MC
85  for (int ii=0; ii<13*2; ii++)
87  }
88 }
double m_minTowerEnergy
Definition: HFClusterAlgo.h:43
double m_maximumRenergy
Definition: HFClusterAlgo.h:43
bool m_forcePulseFlagMC
Definition: HFClusterAlgo.h:45
std::vector< double > m_cutByEta
Definition: HFClusterAlgo.h:48
double m_maximumSL
Definition: HFClusterAlgo.h:43
std::vector< double > m_correctionByEta
Definition: HFClusterAlgo.h:49
std::vector< double > m_seedmnEta
Definition: HFClusterAlgo.h:50
int ii
Definition: cuy.py:588
static const double MCMaterialCorrections_3XX[]
std::vector< double > m_seedMXeta
Definition: HFClusterAlgo.h:51
std::vector< double > m_seedmnPhi
Definition: HFClusterAlgo.h:52
#define M_PI
std::vector< double > m_seedMXphi
Definition: HFClusterAlgo.h:53
double m_seedThreshold
Definition: HFClusterAlgo.h:43

Friends And Related Function Documentation

friend class CompareHFCompleteHitET
friend

Definition at line 40 of file HFClusterAlgo.h.

friend class CompareHFCore
friend

Definition at line 41 of file HFClusterAlgo.h.

Member Data Documentation

std::vector<double> HFClusterAlgo::m_correctionByEta
private

Definition at line 49 of file HFClusterAlgo.h.

int HFClusterAlgo::m_correctionSet
private

Definition at line 47 of file HFClusterAlgo.h.

std::vector<double> HFClusterAlgo::m_cutByEta
private

Definition at line 48 of file HFClusterAlgo.h.

bool HFClusterAlgo::m_forcePulseFlagMC
private

Definition at line 45 of file HFClusterAlgo.h.

bool HFClusterAlgo::m_isMC
private

Definition at line 46 of file HFClusterAlgo.h.

Referenced by isMC().

double HFClusterAlgo::m_maximumRenergy
private

Definition at line 43 of file HFClusterAlgo.h.

double HFClusterAlgo::m_maximumSL
private

Definition at line 43 of file HFClusterAlgo.h.

double HFClusterAlgo::m_minTowerEnergy
private

Definition at line 43 of file HFClusterAlgo.h.

std::vector<double> HFClusterAlgo::m_seedmnEta
private

Definition at line 50 of file HFClusterAlgo.h.

std::vector<double> HFClusterAlgo::m_seedmnPhi
private

Definition at line 52 of file HFClusterAlgo.h.

std::vector<double> HFClusterAlgo::m_seedMXeta
private

Definition at line 51 of file HFClusterAlgo.h.

std::vector<double> HFClusterAlgo::m_seedMXphi
private

Definition at line 53 of file HFClusterAlgo.h.

double HFClusterAlgo::m_seedThreshold
private

Definition at line 43 of file HFClusterAlgo.h.

bool HFClusterAlgo::m_usePMTFlag
private

Definition at line 44 of file HFClusterAlgo.h.

bool HFClusterAlgo::m_usePulseFlag
private

Definition at line 45 of file HFClusterAlgo.h.