CMS 3D CMS Logo

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

#include <CATopJetAlgorithm.h>

Public Member Functions

 CATopJetAlgorithm (edm::InputTag mSrc, bool verbose, int algorithm, int useAdjacency, double centralEtaCut, double ptMin, const std::vector< double > &sumEtBins, const std::vector< double > &rBins, const std::vector< double > &ptFracBins, const std::vector< double > &deltarBins, const std::vector< double > &nCellBins, double seedThreshold, bool useMaxTower, double sumEtEtaCut, double etFrac, boost::shared_ptr< fastjet::JetDefinition > fjJetDefinition, bool doAreaFastjet, boost::shared_ptr< fastjet::ActiveAreaSpec > fjActiveArea, double voronoiRfact)
 
void run (const std::vector< fastjet::PseudoJet > &cell_particles, std::vector< CompoundPseudoJet > &hardjetsOutput)
 Find the ProtoJets from the collection of input Candidates. More...
 

Private Member Functions

bool adjacentCells (const fastjet::PseudoJet &jet1, const fastjet::PseudoJet &jet2, const std::vector< fastjet::PseudoJet > &cell_particles, const fastjet::ClusterSequence &theClusterSequence, double nCellMin) const
 
bool decomposeJet (const fastjet::PseudoJet &theJet, const fastjet::ClusterSequence &theClusterSequence, const std::vector< fastjet::PseudoJet > &cell_particles, double ptHard, double nCellMin, double deltarcut, fastjet::PseudoJet &ja, fastjet::PseudoJet &jb, std::vector< fastjet::PseudoJet > &leftovers) const
 

Private Attributes

int algorithm_
 
double centralEtaCut_
 
std::vector< double > deltarBins_
 
bool doAreaFastjet_
 
double etFrac_
 
boost::shared_ptr
< fastjet::ActiveAreaSpec > 
fjActiveArea_
 
boost::shared_ptr
< fastjet::JetDefinition > 
fjJetDefinition_
 
std::string jetType_
 
edm::InputTag mSrc_
 
std::vector< double > nCellBins_
 
std::vector< double > ptFracBins_
 
double ptMin_
 
std::vector< double > rBins_
 
double seedThreshold_
 
std::vector< double > sumEtBins_
 
double sumEtEtaCut_
 
int useAdjacency_
 
bool useMaxTower_
 
bool verbose_
 
double voronoiRfact_
 

Detailed Description

Definition at line 44 of file CATopJetAlgorithm.h.

Constructor & Destructor Documentation

CATopJetAlgorithm::CATopJetAlgorithm ( edm::InputTag  mSrc,
bool  verbose,
int  algorithm,
int  useAdjacency,
double  centralEtaCut,
double  ptMin,
const std::vector< double > &  sumEtBins,
const std::vector< double > &  rBins,
const std::vector< double > &  ptFracBins,
const std::vector< double > &  deltarBins,
const std::vector< double > &  nCellBins,
double  seedThreshold,
bool  useMaxTower,
double  sumEtEtaCut,
double  etFrac,
boost::shared_ptr< fastjet::JetDefinition >  fjJetDefinition,
bool  doAreaFastjet,
boost::shared_ptr< fastjet::ActiveAreaSpec >  fjActiveArea,
double  voronoiRfact 
)
inline

Constructor

Definition at line 48 of file CATopJetAlgorithm.h.

66  :
67  mSrc_ (mSrc ),
68  verbose_ (verbose ),
70  useAdjacency_ (useAdjacency ),
71  centralEtaCut_ (centralEtaCut ),
72  ptMin_ (ptMin ),
73  sumEtBins_ (sumEtBins ),
74  rBins_ (rBins ),
75  ptFracBins_ (ptFracBins ),
76  deltarBins_ (deltarBins ),
77  nCellBins_ (nCellBins ),
78  seedThreshold_ (seedThreshold ),
79  useMaxTower_ (useMaxTower ),
80  sumEtEtaCut_ (sumEtEtaCut ),
81  etFrac_ (etFrac ),
82  fjJetDefinition_(fjJetDefinition),
83  doAreaFastjet_ (doAreaFastjet),
84  fjActiveArea_ (fjActiveArea),
85  voronoiRfact_ (voronoiRfact)
86  { }
boost::shared_ptr< fastjet::ActiveAreaSpec > fjActiveArea_
< trclass="colgroup">< tdclass="colgroup"colspan=5 > Ecal cluster collections</td ></tr >< tr >< td >< ahref="classreco_1_1BasicCluster.html"> reco::BasicCluster</a ></td >< td >< ahref="DataFormats_EgammaReco.html"> reco::BasicClusterCollection</a ></td >< td >< ahref="#"> hybridSuperClusters</a ></td >< tdclass="description"> Basic clusters reconstructed with hybrid algorithm(barrel only)</td >< td >S.Rahatlou</td ></tr >< tr >< td >< a href
std::vector< double > sumEtBins_
std::vector< double > rBins_
edm::InputTag mSrc_
std::vector< double > ptFracBins_
boost::shared_ptr< fastjet::JetDefinition > fjJetDefinition_
std::vector< double > deltarBins_
std::vector< double > nCellBins_

Member Function Documentation

bool CATopJetAlgorithm::adjacentCells ( const fastjet::PseudoJet &  jet1,
const fastjet::PseudoJet &  jet2,
const std::vector< fastjet::PseudoJet > &  cell_particles,
const fastjet::ClusterSequence &  theClusterSequence,
double  nCellMin 
) const
private

Definition at line 267 of file CATopJetAlgorithm.cc.

References abs, and reco::deltaPhi().

270  {
271 
272 
273  double eta1 = jet1.rapidity();
274  double phi1 = jet1.phi();
275  double eta2 = jet2.rapidity();
276  double phi2 = jet2.phi();
277 
278  double deta = abs(eta2 - eta1) / 0.087;
279  double dphi = fabs( reco::deltaPhi(phi2,phi1) ) / 0.087;
280 
281  return ( ( deta + dphi ) <= nCellMin );
282 }
#define abs(x)
Definition: mlp_lapack.h:159
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
bool CATopJetAlgorithm::decomposeJet ( const fastjet::PseudoJet &  theJet,
const fastjet::ClusterSequence &  theClusterSequence,
const std::vector< fastjet::PseudoJet > &  cell_particles,
double  ptHard,
double  nCellMin,
double  deltarcut,
fastjet::PseudoJet &  ja,
fastjet::PseudoJet &  jb,
std::vector< fastjet::PseudoJet > &  leftovers 
) const
private

Adjacency Requirement ///

Pt Fraction Requirement ///

Definition at line 289 of file CATopJetAlgorithm.cc.

References gather_cfg::cout, Geom::deltaPhi(), deltaR(), and j.

294  {
295 
296  bool goodBreak;
297  fastjet::PseudoJet j = theJet;
298  double InputObjectPt = j.perp();
299  if ( verbose_ )cout<<"Input Object Pt = "<<InputObjectPt<<endl;
300  if ( verbose_ )cout<<"ptHard = "<<ptHard<<endl;
301  leftovers.clear();
302  if ( verbose_ )cout<<"start while loop"<<endl;
303 
304  while (1) { // watch out for infinite loop!
305  goodBreak = theClusterSequence.has_parents(j,ja,jb);
306  if (!goodBreak){
307  if ( verbose_ )cout<<"bad break. this is one cell. can't decluster anymore."<<endl;
308  break; // this is one cell, can't decluster anymore
309  }
310 
311  if ( verbose_ )cout<<"good break. ja Pt = "<<ja.perp()<<" jb Pt = "<<jb.perp()<<endl;
312 
314 
315  // check if clusters are adjacent using a constant deltar adjacency.
316  double clusters_deltar=fabs(ja.eta()-jb.eta())+fabs(deltaPhi(ja.phi(),jb.phi()));
317 
318  if ( verbose_ && useAdjacency_ ==1)cout<<"clusters_deltar = "<<clusters_deltar<<endl;
319  if ( verbose_ && useAdjacency_ ==1)cout<<"deltar cut = "<<deltarcut<<endl;
320 
321  if ( useAdjacency_==1 && clusters_deltar < deltarcut){
322  if ( verbose_ )cout<<"clusters too close. consant adj. break."<<endl;
323  break;
324  }
325 
326  // Check if clusters are adjacent using a DeltaR adjacency which is a function of pT.
327  double clusters_deltaR=deltaR( ja.eta(), ja.phi(), jb.eta(), jb.phi() );
328 
329  if ( verbose_ && useAdjacency_ ==2)cout<<"clusters_deltaR = "<<clusters_deltaR<<endl;
330  if ( verbose_ && useAdjacency_ ==2)cout<<"0.4-0.0004*InputObjectPt = "<<0.4-0.0004*InputObjectPt<<endl;
331 
332  if ( useAdjacency_==2 && clusters_deltaR < 0.4-0.0004*InputObjectPt)
333  {
334  if ( verbose_ )cout<<"clusters too close. modified adj. break."<<endl;
335  break;
336  }
337 
338  // Check if clusters are adjacent in the calorimeter.
339  if ( useAdjacency_==3 && adjacentCells(ja,jb,cell_particles,theClusterSequence,nCellMin) ){
340  if ( verbose_ )cout<<"clusters too close in the calorimeter. calorimeter adj. break."<<endl;
341  break; // the clusters are "adjacent" in the calorimeter => shouldn't have decomposed
342  }
343 
344  if ( verbose_ )cout<<"clusters pass distance cut"<<endl;
345 
347 
348  if ( verbose_ )cout<<"ptHard = "<<ptHard<<endl;
349 
350  if (ja.perp() < ptHard && jb.perp() < ptHard){
351  if ( verbose_ )cout<<"two soft clusters. dead end"<<endl;
352  break; // broke into two soft clusters, dead end
353  }
354 
355  if (ja.perp() > ptHard && jb.perp() > ptHard){
356  if ( verbose_ )cout<<"two hard clusters. done"<<endl;
357  return true; // broke into two hard clusters, we're done!
358  }
359 
360  else if (ja.perp() > jb.perp()) { // broke into one hard and one soft, ditch the soft one and try again
361  if ( verbose_ )cout<<"ja hard jb soft. try to split hard. j = ja"<<endl;
362  j = ja;
363  vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(jb);
364  leftovers.insert(leftovers.end(),particles.begin(),particles.end());
365  }
366  else {
367  if ( verbose_ )cout<<"ja hard jb soft. try to split hard. j = jb"<<endl;
368  j = jb;
369  vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(ja);
370  leftovers.insert(leftovers.end(),particles.begin(),particles.end());
371  }
372  }
373 
374  if ( verbose_ )cout<<"did not decluster."<<endl; // did not decluster into hard subjets
375 
376  ja.reset(0,0,0,0);
377  jb.reset(0,0,0,0);
378  leftovers.clear();
379  return false;
380 }
double deltaPhi(float phi1, float phi2)
Definition: VectorUtil.h:30
bool adjacentCells(const fastjet::PseudoJet &jet1, const fastjet::PseudoJet &jet2, const std::vector< fastjet::PseudoJet > &cell_particles, const fastjet::ClusterSequence &theClusterSequence, double nCellMin) const
int j
Definition: DBlmapReader.cc:9
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
tuple cout
Definition: gather_cfg.py:41
void CATopJetAlgorithm::run ( const std::vector< fastjet::PseudoJet > &  cell_particles,
std::vector< CompoundPseudoJet > &  hardjetsOutput 
)

Find the ProtoJets from the collection of input Candidates.

Definition at line 20 of file CATopJetAlgorithm.cc.

References gather_cfg::cout, i, perp(), and python.multivaluedict::sort().

Referenced by cms::CATopJetProducer::runAlgorithm().

22 {
23  if ( verbose_ ) cout << "Welcome to CATopSubJetAlgorithm::run" << endl;
24 
25  // Sum Et of the event
26  double sumEt = 0.;
27 
28  //make a list of input objects ordered by ET and calculate sum et
29  // list of fastjet pseudojet constituents
30  for (unsigned i = 0; i < cell_particles.size(); ++i) {
31  sumEt += cell_particles[i].perp();
32  }
33 
34  // Determine which bin we are in for et clustering
35  int sumEtBinId = -1;
36  for ( unsigned int i = 0; i < sumEtBins_.size(); ++i ) {
37  if ( sumEt > sumEtBins_[i] ) sumEtBinId = i;
38  }
39  if ( verbose_ ) cout << "Using sumEt = " << sumEt << ", bin = " << sumEtBinId << endl;
40 
41  // If the sum et is too low, exit
42  if ( sumEtBinId < 0 ) {
43  return;
44  }
45 
46  // empty 4-vector
47  fastjet::PseudoJet blankJetA(0,0,0,0);
48  blankJetA.set_user_index(-1);
49  const fastjet::PseudoJet blankJet = blankJetA;
50 
51  // Define adjacency variables which depend on which sumEtBin we are in
52  double deltarcut = deltarBins_[sumEtBinId];
53  double nCellMin = nCellBins_[sumEtBinId];
54 
55  if ( verbose_ )cout<<"useAdjacency_ = "<<useAdjacency_<<endl;
56  if ( verbose_ && useAdjacency_==0)cout<<"No Adjacency"<<endl;
57  if ( verbose_ && useAdjacency_==1)cout<<"using deltar adjacency"<<endl;
58  if ( verbose_ && useAdjacency_==2)cout<<"using modified adjacency"<<endl;
59  if ( verbose_ && useAdjacency_==3)cout<<"using calorimeter nearest neighbor based adjacency"<<endl;
60  if ( verbose_ && useAdjacency_==1)cout << "Using deltarcut = " << deltarcut << endl;
61  if ( verbose_ && useAdjacency_==3)cout << "Using nCellMin = " << nCellMin << endl;
62 
63 
64  // Define strategy, recombination scheme, and jet definition
65  fastjet::JetDefinition jetDef( fjJetDefinition_->jet_algorithm(),
66  rBins_[sumEtBinId],
67  fjJetDefinition_->recombination_scheme(),
68  fjJetDefinition_->strategy() );
69 
70  if ( verbose_ ) cout << "About to do jet clustering in CA" << endl;
71  // run the jet clustering
72 
73  //cluster the jets with the jet definition jetDef:
74  // run algorithm
75  boost::shared_ptr<fastjet::ClusterSequence> fjClusterSeq;
76  if ( !doAreaFastjet_ ) {
77  fjClusterSeq = boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequence( cell_particles, jetDef ) );
78  } else if (voronoiRfact_ <= 0) {
79  fjClusterSeq = boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequenceArea( cell_particles, jetDef , *fjActiveArea_ ) );
80  } else {
81  fjClusterSeq = boost::shared_ptr<fastjet::ClusterSequence>( new fastjet::ClusterSequenceVoronoiArea( cell_particles, jetDef , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
82  }
83 
84  if ( verbose_ ) cout << "Getting inclusive jets" << endl;
85  // Get the transient inclusive jets
86  vector<fastjet::PseudoJet> inclusiveJets = fjClusterSeq->inclusive_jets(ptMin_);
87 
88  if ( verbose_ ) cout << "Getting central jets" << endl;
89  // Find the transient central jets
90  vector<fastjet::PseudoJet> centralJets;
91  for (unsigned int i = 0; i < inclusiveJets.size(); i++) {
92 
93  if (inclusiveJets[i].perp() > ptMin_ && fabs(inclusiveJets[i].rapidity()) < centralEtaCut_) {
94  centralJets.push_back(inclusiveJets[i]);
95  }
96  }
97  // Sort the transient central jets in Et
98  GreaterByEtPseudoJet compEt;
99  sort( centralJets.begin(), centralJets.end(), compEt );
100 
101  // These will store the 4-vectors of each hard jet
102  vector<math::XYZTLorentzVector> p4_hardJets;
103 
104  // These will store the indices of each subjet that
105  // are present in each jet
106  vector<vector<int> > indices( centralJets.size() );
107 
108  // Loop over central jets, attempt to find substructure
109  vector<fastjet::PseudoJet>::iterator jetIt = centralJets.begin(),
110  centralJetsBegin = centralJets.begin(),
111  centralJetsEnd = centralJets.end();
112  if ( verbose_ )cout<<"Loop over jets"<<endl;
113  int i=0;
114  for ( ; jetIt != centralJetsEnd; ++jetIt ) {
115  if ( verbose_ )cout<<"\nJet "<<i<<endl;
116  i++;
117  fastjet::PseudoJet localJet = *jetIt;
118 
119  // Get the 4-vector for this jet
120  p4_hardJets.push_back( math::XYZTLorentzVector(localJet.px(), localJet.py(), localJet.pz(), localJet.e() ));
121 
122  // jet decomposition. try to find 3 or 4 hard, well-localized subjets, characteristic of a boosted top.
123  if ( verbose_ )cout<<"local jet pt = "<<localJet.perp()<<endl;
124  if ( verbose_ )cout<<"deltap = "<<ptFracBins_[sumEtBinId]<<endl;
125 
126  double ptHard = ptFracBins_[sumEtBinId]*localJet.perp();
127  vector<fastjet::PseudoJet> leftoversAll;
128 
129  // stage 1: primary decomposition. look for when the jet declusters into two hard subjets
130  if ( verbose_ ) cout << "Doing decomposition 1" << endl;
131  fastjet::PseudoJet ja, jb;
132  vector<fastjet::PseudoJet> leftovers1;
133  bool hardBreak1 = decomposeJet(localJet,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,ja,jb,leftovers1);
134  leftoversAll.insert(leftoversAll.end(),leftovers1.begin(),leftovers1.end());
135 
136  // stage 2: secondary decomposition. look for when the hard subjets found above further decluster into two hard sub-subjets
137  //
138  // ja -> jaa+jab ?
139  if ( verbose_ ) cout << "Doing decomposition 2. ja->jaa+jab?" << endl;
140  fastjet::PseudoJet jaa, jab;
141  vector<fastjet::PseudoJet> leftovers2a;
142  bool hardBreak2a = false;
143  if (hardBreak1) hardBreak2a = decomposeJet(ja,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jaa,jab,leftovers2a);
144  leftoversAll.insert(leftoversAll.end(),leftovers2a.begin(),leftovers2a.end());
145  // jb -> jba+jbb ?
146  if ( verbose_ ) cout << "Doing decomposition 2. ja->jba+jbb?" << endl;
147  fastjet::PseudoJet jba, jbb;
148  vector<fastjet::PseudoJet> leftovers2b;
149  bool hardBreak2b = false;
150  if (hardBreak1) hardBreak2b = decomposeJet(jb,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jba,jbb,leftovers2b);
151  leftoversAll.insert(leftoversAll.end(),leftovers2b.begin(),leftovers2b.end());
152 
153  // NOTE: it might be good to consider some checks for whether these subjets can be further decomposed. e.g., the above procedure leaves
154  // open the possibility of "subjets" that actually consist of two or more distinct hard clusters. however, this kind of thing
155  // is a rarity for the simulations so far considered.
156 
157  // proceed if one or both of the above hard subjets successfully decomposed
158  if ( verbose_ ) cout << "Done with decomposition" << endl;
159 
160  if ( verbose_ ) cout<<"hardBreak1 = "<<hardBreak1<<endl;
161  if ( verbose_ ) cout<<"hardBreak2a = "<<hardBreak2a<<endl;
162  if ( verbose_ ) cout<<"hardBreak2b = "<<hardBreak2b<<endl;
163 
164  fastjet::PseudoJet hardA = blankJet, hardB = blankJet, hardC = blankJet, hardD = blankJet;
165  if (!hardBreak1) {
166  hardA = localJet;
167  hardB = blankJet;
168  hardC = blankJet;
169  hardD = blankJet;
170  if(verbose_)cout<<"Hardbreak failed. Save subjet1=localJet"<<endl;
171  }
172  if (hardBreak1 && !hardBreak2a && !hardBreak2b) {
173  hardA = ja;
174  hardB = jb;
175  hardC = blankJet;
176  hardD = blankJet;
177  if(verbose_)cout<<"First decomposition succeeded, both second decompositions failed. Save subjet1=ja subjet2=jb"<<endl;
178  }
179  if (hardBreak1 && hardBreak2a && !hardBreak2b) {
180  hardA = jaa;
181  hardB = jab;
182  hardC = jb;
183  hardD = blankJet;
184  if(verbose_)cout<<"First decomposition succeeded, ja split succesfully, jb did not split. Save subjet1=jaa subjet2=jab subjet3=jb"<<endl;
185  }
186  if (hardBreak1 && !hardBreak2a && hardBreak2b) {
187  hardA = jba;
188  hardB = jbb;
189  hardC = ja;
190  hardD = blankJet;
191  if(verbose_)cout<<"First decomposition succeeded, jb split succesfully, ja did not split. Save subjet1=jba subjet2=jbb subjet3=ja"<<endl;
192  }
193  if (hardBreak1 && hardBreak2a && hardBreak2b) {
194  hardA = jaa;
195  hardB = jab;
196  hardC = jba;
197  hardD = jbb;
198  if(verbose_)cout<<"First decomposition and both secondary decompositions succeeded. Save subjet1=jaa subjet2=jab subjet3=jba subjet4=jbb"<<endl;
199  }
200 
201  // check if we are left with >= 3 hard subjets
202  fastjet::PseudoJet subjet1 = blankJet;
203  fastjet::PseudoJet subjet2 = blankJet;
204  fastjet::PseudoJet subjet3 = blankJet;
205  fastjet::PseudoJet subjet4 = blankJet;
206  subjet1 = hardA; subjet2 = hardB; subjet3 = hardC; subjet4 = hardD;
207 
208  // record the hard subjets
209  vector<fastjet::PseudoJet> hardSubjets;
210 
211  if ( subjet1.user_index() >= 0 )
212  hardSubjets.push_back(subjet1);
213  if ( subjet2.user_index() >= 0 )
214  hardSubjets.push_back(subjet2);
215  if ( subjet3.user_index() >= 0 )
216  hardSubjets.push_back(subjet3);
217  if ( subjet4.user_index() >= 0 )
218  hardSubjets.push_back(subjet4);
219  sort(hardSubjets.begin(), hardSubjets.end(), compEt );
220 
221  // create the subjets objects to put into the "output" objects
222  vector<CompoundPseudoSubJet> subjetsOutput;
223  std::vector<fastjet::PseudoJet>::const_iterator itSubJetBegin = hardSubjets.begin(),
224  itSubJet = itSubJetBegin, itSubJetEnd = hardSubjets.end();
225  for (; itSubJet != itSubJetEnd; ++itSubJet ){
226  // if ( verbose_ ) cout << "Adding input collection element " << (*itSubJet).user_index() << endl;
227  // if ( (*itSubJet).user_index() >= 0 && (*itSubJet).user_index() < cell_particles.size() )
228 
229  // Get the transient subjet constituents from fastjet
230  vector<fastjet::PseudoJet> subjetFastjetConstituents = fjClusterSeq->constituents( *itSubJet );
231 
232  // Get the indices of the constituents
233  vector<int> constituents;
234 
235  // Loop over the constituents and get the indices
236  vector<fastjet::PseudoJet>::const_iterator fastSubIt = subjetFastjetConstituents.begin(),
237  transConstBegin = subjetFastjetConstituents.begin(),
238  transConstEnd = subjetFastjetConstituents.end();
239  for ( ; fastSubIt != transConstEnd; ++fastSubIt ) {
240  if ( fastSubIt->user_index() >= 0 && static_cast<unsigned int>(fastSubIt->user_index()) < cell_particles.size() ) {
241  constituents.push_back( fastSubIt->user_index() );
242  }
243  }
244 
245  // Make a CompoundPseudoSubJet object to hold this subjet and the indices of its constituents
246  subjetsOutput.push_back( CompoundPseudoSubJet( *itSubJet, constituents ) );
247  }
248 
249 
250  double fatJetArea = (doAreaFastjet_) ?
251  ((fastjet::ClusterSequenceArea&)*fjClusterSeq).area(*jetIt) : 0.0;
252 
253  // Make a CompoundPseudoJet object to hold this hard jet, and the subjets that make it up
254  hardjetsOutput.push_back( CompoundPseudoJet( *jetIt,fatJetArea,subjetsOutput));
255  }
256 }
boost::shared_ptr< fastjet::ActiveAreaSpec > fjActiveArea_
int i
Definition: DBlmapReader.cc:9
std::vector< double > sumEtBins_
std::vector< double > rBins_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
bool decomposeJet(const fastjet::PseudoJet &theJet, const fastjet::ClusterSequence &theClusterSequence, const std::vector< fastjet::PseudoJet > &cell_particles, double ptHard, double nCellMin, double deltarcut, fastjet::PseudoJet &ja, fastjet::PseudoJet &jb, std::vector< fastjet::PseudoJet > &leftovers) const
CompoundPseudoJet holds an association of fastjet::PseudoJets that represent a &quot;hard&quot; top jet with su...
std::vector< double > ptFracBins_
boost::shared_ptr< fastjet::JetDefinition > fjJetDefinition_
std::vector< double > deltarBins_
T perp() const
Magnitude of transverse component.
std::vector< double > nCellBins_
tuple cout
Definition: gather_cfg.py:41

Member Data Documentation

int CATopJetAlgorithm::algorithm_
private

Definition at line 97 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::centralEtaCut_
private

Definition at line 103 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::deltarBins_
private

Definition at line 108 of file CATopJetAlgorithm.h.

bool CATopJetAlgorithm::doAreaFastjet_
private

Definition at line 117 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::etFrac_
private

Definition at line 114 of file CATopJetAlgorithm.h.

boost::shared_ptr<fastjet::ActiveAreaSpec> CATopJetAlgorithm::fjActiveArea_
private

Definition at line 118 of file CATopJetAlgorithm.h.

boost::shared_ptr<fastjet::JetDefinition> CATopJetAlgorithm::fjJetDefinition_
private

Definition at line 116 of file CATopJetAlgorithm.h.

std::string CATopJetAlgorithm::jetType_
private

Definition at line 115 of file CATopJetAlgorithm.h.

edm::InputTag CATopJetAlgorithm::mSrc_
private

Definition at line 95 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::nCellBins_
private

Definition at line 109 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::ptFracBins_
private

Definition at line 107 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::ptMin_
private

Definition at line 104 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::rBins_
private

Definition at line 106 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::seedThreshold_
private

Definition at line 111 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::sumEtBins_
private

Definition at line 105 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::sumEtEtaCut_
private

Definition at line 113 of file CATopJetAlgorithm.h.

int CATopJetAlgorithm::useAdjacency_
private

Definition at line 98 of file CATopJetAlgorithm.h.

bool CATopJetAlgorithm::useMaxTower_
private

Definition at line 112 of file CATopJetAlgorithm.h.

bool CATopJetAlgorithm::verbose_
private

Definition at line 96 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::voronoiRfact_
private

Definition at line 119 of file CATopJetAlgorithm.h.