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)
 
void run (const std::vector< fastjet::PseudoJet > &cell_particles, std::vector< fastjet::PseudoJet > &hardjetsOutput, boost::shared_ptr< fastjet::ClusterSequence > &fjClusterSeq)
 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_
 
double etFrac_
 
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_
 

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 
)
inline

Constructor

Definition at line 48 of file CATopJetAlgorithm.h.

62  :
63  mSrc_ (mSrc ),
64  verbose_ (verbose ),
66  useAdjacency_ (useAdjacency ),
67  centralEtaCut_ (centralEtaCut ),
68  ptMin_ (ptMin ),
69  sumEtBins_ (sumEtBins ),
70  rBins_ (rBins ),
71  ptFracBins_ (ptFracBins ),
72  deltarBins_ (deltarBins ),
73  nCellBins_ (nCellBins ),
74  seedThreshold_ (seedThreshold ),
75  useMaxTower_ (useMaxTower ),
76  sumEtEtaCut_ (sumEtEtaCut ),
77  etFrac_ (etFrac )
78 
79  { }
< 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_
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 269 of file CATopJetAlgorithm.cc.

References abs, and reco::deltaPhi().

272  {
273 
274 
275  double eta1 = jet1.rapidity();
276  double phi1 = jet1.phi();
277  double eta2 = jet2.rapidity();
278  double phi2 = jet2.phi();
279 
280  double deta = abs(eta2 - eta1) / 0.087;
281  double dphi = fabs( reco::deltaPhi(phi2,phi1) ) / 0.087;
282 
283  return ( ( deta + dphi ) <= nCellMin );
284 }
#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 291 of file CATopJetAlgorithm.cc.

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

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

Find the ProtoJets from the collection of input Candidates.

Definition at line 20 of file CATopJetAlgorithm.cc.

References gather_cfg::cout, i, errorMatrix2Lands::indices, join(), perp(), csv2json::pieces, and python.multivaluedict::sort().

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

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

Member Data Documentation

int CATopJetAlgorithm::algorithm_
private

Definition at line 91 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::centralEtaCut_
private

Definition at line 97 of file CATopJetAlgorithm.h.

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

Definition at line 102 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::etFrac_
private

Definition at line 108 of file CATopJetAlgorithm.h.

std::string CATopJetAlgorithm::jetType_
private

Definition at line 109 of file CATopJetAlgorithm.h.

edm::InputTag CATopJetAlgorithm::mSrc_
private

Definition at line 89 of file CATopJetAlgorithm.h.

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

Definition at line 103 of file CATopJetAlgorithm.h.

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

Definition at line 101 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::ptMin_
private

Definition at line 98 of file CATopJetAlgorithm.h.

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

Definition at line 100 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::seedThreshold_
private

Definition at line 105 of file CATopJetAlgorithm.h.

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

Definition at line 99 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::sumEtEtaCut_
private

Definition at line 107 of file CATopJetAlgorithm.h.

int CATopJetAlgorithm::useAdjacency_
private

Definition at line 92 of file CATopJetAlgorithm.h.

bool CATopJetAlgorithm::useMaxTower_
private

Definition at line 106 of file CATopJetAlgorithm.h.

bool CATopJetAlgorithm::verbose_
private

Definition at line 90 of file CATopJetAlgorithm.h.