CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastjetJetProducer.cc
Go to the documentation of this file.
1 //
3 // FastjetJetProducer
4 // ------------------
5 //
6 // 04/21/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
8 
10 
19 
25 
28 
29 #include "fastjet/SISConePlugin.hh"
30 #include "fastjet/CMSIterativeConePlugin.hh"
31 #include "fastjet/ATLASConePlugin.hh"
32 #include "fastjet/CDFMidPointPlugin.hh"
33 #include "fastjet/tools/Filter.hh"
34 #include "fastjet/tools/Pruner.hh"
35 #include "fastjet/tools/MassDropTagger.hh"
36 #include "fastjet/contrib/SoftDrop.hh"
37 #include "fastjet/tools/JetMedianBackgroundEstimator.hh"
38 #include "fastjet/contrib/ConstituentSubtractor.hh"
40 
41 #include <iostream>
42 #include <memory>
43 #include <algorithm>
44 #include <limits>
45 #include <cmath>
46 //#include <fstream>
47 
48 using namespace std;
49 
50 
51 
53 // construction / destruction
55 
56 //______________________________________________________________________________
58  : VirtualJetProducer( iConfig ),
59  useMassDropTagger_(false),
60  useFiltering_(false),
61  useDynamicFiltering_(false),
62  useTrimming_(false),
63  usePruning_(false),
64  useCMSBoostedTauSeedingAlgorithm_(false),
65  useKtPruning_(false),
66  useConstituentSubtraction_(false),
67  useSoftDrop_(false),
68  muCut_(-1.0),
69  yCut_(-1.0),
70  rFilt_(-1.0),
71  rFiltFactor_(-1.0),
72  nFilt_(-1),
73  trimPtFracMin_(-1.0),
74  zCut_(-1.0),
75  RcutFactor_(-1.0),
76  csRho_EtaMax_(-1.0),
77  csRParam_(-1.0),
78  beta_(-1.0)
79 {
80 
81  if ( iConfig.exists("UseOnlyVertexTracks") )
82  useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
83  else
84  useOnlyVertexTracks_ = false;
85 
86  if ( iConfig.exists("UseOnlyOnePV") )
87  useOnlyOnePV_ = iConfig.getParameter<bool>("UseOnlyOnePV");
88  else
89  useOnlyOnePV_ = false;
90 
91  if ( iConfig.exists("DzTrVtxMax") )
92  dzTrVtxMax_ = iConfig.getParameter<double>("DzTrVtxMax");
93  else
94  dzTrVtxMax_ = 999999.;
95  if ( iConfig.exists("DxyTrVtxMax") )
96  dxyTrVtxMax_ = iConfig.getParameter<double>("DxyTrVtxMax");
97  else
98  dxyTrVtxMax_ = 999999.;
99  if ( iConfig.exists("MinVtxNdof") )
100  minVtxNdof_ = iConfig.getParameter<int>("MinVtxNdof");
101  else
102  minVtxNdof_ = 5;
103  if ( iConfig.exists("MaxVtxZ") )
104  maxVtxZ_ = iConfig.getParameter<double>("MaxVtxZ");
105  else
106  maxVtxZ_ = 15;
107 
108 
109  if ( iConfig.exists("useFiltering") ||
110  iConfig.exists("useTrimming") ||
111  iConfig.exists("usePruning") ||
112  iConfig.exists("useMassDropTagger") ||
113  iConfig.exists("useCMSBoostedTauSeedingAlgorithm") ||
114  iConfig.exists("useConstituentSubtraction") ||
115  iConfig.exists("useSoftDrop")
116  ) {
117  useMassDropTagger_=false;
118  useFiltering_=false;
119  useDynamicFiltering_=false;
120  useTrimming_=false;
121  usePruning_=false;
123  useKtPruning_=false;
125  useSoftDrop_ = false;
126  rFilt_=-1.0;
127  rFiltFactor_=-1.0;
128  nFilt_=-1;
129  trimPtFracMin_=-1.0;
130  zCut_=-1.0;
131  RcutFactor_=-1.0;
132  muCut_=-1.0;
133  yCut_=-1.0;
134  subjetPtMin_ = -1.0;
135  muMin_ = -1.0;
136  muMax_ = -1.0;
137  yMin_ = -1.0;
138  yMax_ = -1.0;
139  dRMin_ = -1.0;
140  dRMax_ = -1.0;
141  maxDepth_ = -1;
142  csRho_EtaMax_ = -1.0;
143  csRParam_ = -1.0;
144  beta_ = -1.0;
145  useExplicitGhosts_ = true;
146 
147  if ( iConfig.exists("useMassDropTagger") ) {
148  useMassDropTagger_ = iConfig.getParameter<bool>("useMassDropTagger");
149  muCut_ = iConfig.getParameter<double>("muCut");
150  yCut_ = iConfig.getParameter<double>("yCut");
151  }
152 
153  if ( iConfig.exists("useFiltering") ) {
154  useFiltering_ = iConfig.getParameter<bool>("useFiltering");
155  rFilt_ = iConfig.getParameter<double>("rFilt");
156  nFilt_ = iConfig.getParameter<int>("nFilt");
157  if ( iConfig.exists("useDynamicFiltering") ) {
158  useDynamicFiltering_ = iConfig.getParameter<bool>("useDynamicFiltering");
159  rFiltFactor_ = iConfig.getParameter<double>("rFiltFactor");
160  if ( useDynamicFiltering_ )
161  rFiltDynamic_ = DynamicRfiltPtr(new DynamicRfilt(rFilt_, rFiltFactor_));
162  }
163  }
164 
165  if ( iConfig.exists("useTrimming") ) {
166  useTrimming_ = iConfig.getParameter<bool>("useTrimming");
167  rFilt_ = iConfig.getParameter<double>("rFilt");
168  trimPtFracMin_ = iConfig.getParameter<double>("trimPtFracMin");
169  }
170 
171  if ( iConfig.exists("usePruning") ) {
172  usePruning_ = iConfig.getParameter<bool>("usePruning");
173  zCut_ = iConfig.getParameter<double>("zcut");
174  RcutFactor_ = iConfig.getParameter<double>("rcut_factor");
175  nFilt_ = iConfig.getParameter<int>("nFilt");
176  if ( iConfig.exists("useKtPruning") )
177  useKtPruning_ = iConfig.getParameter<bool>("useKtPruning");
178  }
179 
180  if ( iConfig.exists("useCMSBoostedTauSeedingAlgorithm") ) {
181  useCMSBoostedTauSeedingAlgorithm_ = iConfig.getParameter<bool>("useCMSBoostedTauSeedingAlgorithm");
182  subjetPtMin_ = iConfig.getParameter<double>("subjetPtMin");
183  muMin_ = iConfig.getParameter<double>("muMin");
184  muMax_ = iConfig.getParameter<double>("muMax");
185  yMin_ = iConfig.getParameter<double>("yMin");
186  yMax_ = iConfig.getParameter<double>("yMax");
187  dRMin_ = iConfig.getParameter<double>("dRMin");
188  dRMax_ = iConfig.getParameter<double>("dRMax");
189  maxDepth_ = iConfig.getParameter<int>("maxDepth");
190  }
191 
192  if ( iConfig.exists("useConstituentSubtraction") ) {
193 
194  if ( fjAreaDefinition_.get() == 0 ) {
195  throw cms::Exception("AreaMustBeSet") << "Logic error. The area definition must be set if you use constituent subtraction." << std::endl;
196  }
197 
198  useConstituentSubtraction_ = iConfig.getParameter<bool>("useConstituentSubtraction");
199  csRho_EtaMax_ = iConfig.getParameter<double>("csRho_EtaMax");
200  csRParam_ = iConfig.getParameter<double>("csRParam");
201  }
202 
203  if ( iConfig.exists("useSoftDrop") ) {
204  if ( usePruning_ ) {
205  throw cms::Exception("PruningAndSoftDrop") << "Logic error. Soft drop is a generalized pruning, do not run them together." << std::endl;
206  }
207  useSoftDrop_ = iConfig.getParameter<bool>("useSoftDrop");
208  zCut_ = iConfig.getParameter<double>("zcut");
209  beta_ = iConfig.getParameter<double>("beta");
210  }
211 
212  }
213 
214  input_chrefcand_token_ = consumes<edm::View<reco::RecoChargedRefCandidate> >(src_);
215 
216 }
217 
218 
219 //______________________________________________________________________________
221 {
222 }
223 
224 
226 // implementation of member functions
228 
230 {
231 
232  // for everything but track jets
233  if (!makeTrackJet(jetTypeE)) {
234 
235  // use the default production from one collection
236  VirtualJetProducer::produce( iEvent, iSetup );
237 
238  } else { // produce trackjets from tracks grouped per primary vertex
239 
240  produceTrackJets(iEvent, iSetup);
241 
242  }
243 
244 }
245 
246 
248 {
249 
250  // read in the track candidates
252  iEvent.getByToken(input_chrefcand_token_, inputsHandle);
253 
254  // make collection with pointers so we can play around with it
255  std::vector<edm::Ptr<reco::RecoChargedRefCandidate> > allInputs;
256  std::vector<edm::Ptr<reco::Candidate> > origInputs;
257  for (size_t i = 0; i < inputsHandle->size(); ++i) {
258  allInputs.push_back(inputsHandle->ptrAt(i));
259  origInputs.push_back(inputsHandle->ptrAt(i));
260  }
261 
262  // read in the PV collection
264  iEvent.getByToken(input_vertex_token_, pvCollection);
265  // define the overall output jet container
266  std::auto_ptr<std::vector<reco::TrackJet> > jets(new std::vector<reco::TrackJet>() );
267 
268  // loop over the good vertices, clustering for each vertex separately
269  for (reco::VertexCollection::const_iterator itVtx = pvCollection->begin(); itVtx != pvCollection->end(); ++itVtx) {
270  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
271 
272  // clear the intermediate containers
273  inputs_.clear();
274  fjInputs_.clear();
275  fjJets_.clear();
276 
277  // if only vertex-associated tracks should be used
278  if (useOnlyVertexTracks_) {
279  // loop over the tracks associated to the vertex
280  for (reco::Vertex::trackRef_iterator itTr = itVtx->tracks_begin(); itTr != itVtx->tracks_end(); ++itTr) {
281  // whether a match was found in the track candidate input
282  bool found = false;
283  // loop over input track candidates
284  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
285  // match the input track candidate to the track from the vertex
286  reco::TrackRef trref(itTr->castTo<reco::TrackRef>());
287  // check if the tracks match
288  if ((*itIn)->track() == trref) {
289  found = true;
290  // add this track candidate to the input for clustering
291  inputs_.push_back(*itIn);
292  // erase the track candidate from the total list of input, so we don't reuse it later
293  allInputs.erase(itIn);
294  // found the candidate track corresponding to the vertex track, so stop the loop
295  break;
296  } // end if match found
297  } // end loop over input tracks
298  // give an info message in case no match is found (can happen if candidates are subset of tracks used for clustering)
299  if (!found) edm::LogInfo("FastjetTrackJetProducer") << "Ignoring a track at vertex which is not in input track collection!";
300  } // end loop over tracks associated to vertex
301  // if all inpt track candidates should be used
302  } else {
303  // loop over input track candidates
304  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
305  // check if the track is close enough to the vertex
306  float dz = (*itIn)->track()->dz(itVtx->position());
307  float dxy = (*itIn)->track()->dxy(itVtx->position());
308  if (fabs(dz) > dzTrVtxMax_) continue;
309  if (fabs(dxy) > dxyTrVtxMax_) continue;
310  bool closervtx = false;
311  // now loop over the good vertices a second time
312  for (reco::VertexCollection::const_iterator itVtx2 = pvCollection->begin(); itVtx2 != pvCollection->end(); ++itVtx2) {
313  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
314  // and check this track is closer to any other vertex (if more than 1 vertex considered)
315  if (!useOnlyOnePV_ &&
316  itVtx != itVtx2 &&
317  fabs((*itIn)->track()->dz(itVtx2->position())) < fabs(dz)) {
318  closervtx = true;
319  break; // 1 closer vertex makes the track already not matched, so break
320  }
321  }
322  // don't add this track if another vertex is found closer
323  if (closervtx) continue;
324  // add this track candidate to the input for clustering
325  inputs_.push_back(*itIn);
326  // erase the track candidate from the total list of input, so we don't reuse it later
327  allInputs.erase(itIn);
328  // take a step back in the loop since we just erased
329  --itIn;
330  }
331  }
332 
333  // convert candidates in inputs_ to fastjet::PseudoJets in fjInputs_
334  fjInputs_.reserve(inputs_.size());
335  inputTowers();
336  LogDebug("FastjetTrackJetProducer") << "Inputted towers\n";
337 
338  // run algorithm, using fjInputs_, modifying fjJets_ and allocating fjClusterSeq_
339  runAlgorithm(iEvent, iSetup);
340  LogDebug("FastjetTrackJetProducer") << "Ran algorithm\n";
341 
342  // convert our jets and add to the overall jet vector
343  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
344  // get the constituents from fastjet
345  std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJets_[ijet]));
346  // convert them to CandidatePtr vector
347  std::vector<reco::CandidatePtr> constituents = getConstituents(fjConstituents);
348  // fill the trackjet
350  // write the specifics to the jet (simultaneously sets 4-vector, vertex).
351  writeSpecific( jet,
352  reco::Particle::LorentzVector(fjJets_[ijet].px(), fjJets_[ijet].py(), fjJets_[ijet].pz(), fjJets_[ijet].E()),
353  vertex_, constituents, iSetup);
354  jet.setJetArea(0);
355  jet.setPileup(0);
356  jet.setPrimaryVertex(edm::Ref<reco::VertexCollection>(pvCollection, (int) (itVtx-pvCollection->begin())));
357  jet.setVertex(itVtx->position());
358  jets->push_back(jet);
359  }
360 
361  if (useOnlyOnePV_) break; // stop vertex loop if only one vertex asked for
362  } // end loop over vertices
363 
364  // put the jets in the collection
365  LogDebug("FastjetTrackJetProducer") << "Put " << jets->size() << " jets in the event.\n";
366  iEvent.put(jets);
367 
368 }
369 
370 
371 //______________________________________________________________________________
373 {
374  // run algorithm
375  /*
376  fjInputs_.clear();
377  double px, py , pz, E;
378  string line;
379  std::ifstream fin("dump3.txt");
380  while (getline(fin, line)){
381  if (line == "#END") break;
382  if (line.substr(0,1) == "#") {continue;}
383  istringstream istr(line);
384  istr >> px >> py >> pz >> E;
385  // create a fastjet::PseudoJet with these components and put it onto
386  // back of the input_particles vector
387  fastjet::PseudoJet j(px,py,pz,E);
388  //if ( fabs(j.rap()) < inputEtaMax )
389  fjInputs_.push_back(fastjet::PseudoJet(px,py,pz,E));
390  }
391  fin.close();
392  */
393 
394  if ( !doAreaFastjet_ && !doRhoFastjet_) {
395  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
396  } else if (voronoiRfact_ <= 0) {
397  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjAreaDefinition_ ) );
398  } else {
399  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
400  }
401 
403  fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
404  } else {
405  fjJets_.clear();
406 
407 
408  transformer_coll transformers;
409 
410 
411  std::vector<fastjet::PseudoJet> tempJets = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
412 
413  unique_ptr<fastjet::JetMedianBackgroundEstimator> bge_rho;
415  fastjet::Selector rho_range = fastjet::SelectorAbsRapMax(csRho_EtaMax_);
416  bge_rho = unique_ptr<fastjet::JetMedianBackgroundEstimator> (new fastjet::JetMedianBackgroundEstimator(rho_range, fastjet::JetDefinition(fastjet::kt_algorithm, csRParam_), *fjAreaDefinition_) );
417  bge_rho->set_particles(fjInputs_);
418  fastjet::contrib::ConstituentSubtractor * constituentSubtractor = new fastjet::contrib::ConstituentSubtractor(bge_rho.get());
419  // this sets the same background estimator to be used for deltaMass density, rho_m, as for pt density, rho:
420  constituentSubtractor->use_common_bge_for_rho_and_rhom(true); // for massless input particles it does not make any difference (rho_m is always zero)
421 
422  transformers.push_back( transformer_ptr(constituentSubtractor) );
423  };
424  if ( useMassDropTagger_ ) {
425  fastjet::MassDropTagger * md_tagger = new fastjet::MassDropTagger ( muCut_, yCut_ );
426  transformers.push_back( transformer_ptr(md_tagger) );
427  }
429  fastjet::contrib::CMSBoostedTauSeedingAlgorithm * tau_tagger =
430  new fastjet::contrib::CMSBoostedTauSeedingAlgorithm ( subjetPtMin_, muMin_, muMax_, yMin_, yMax_, dRMin_, dRMax_, maxDepth_, verbosity_ );
431  transformers.push_back( transformer_ptr(tau_tagger ));
432  }
433  if ( useTrimming_ ) {
434  fastjet::Filter * trimmer = new fastjet::Filter(fastjet::JetDefinition(fastjet::kt_algorithm, rFilt_), fastjet::SelectorPtFractionMin(trimPtFracMin_));
435  transformers.push_back( transformer_ptr(trimmer) );
436  }
437  if ( (useFiltering_) && (!useDynamicFiltering_) ) {
438  fastjet::Filter * filter = new fastjet::Filter(fastjet::JetDefinition(fastjet::cambridge_algorithm, rFilt_), fastjet::SelectorNHardest(nFilt_));
439  transformers.push_back( transformer_ptr(filter));
440  }
441 
442  if ( (usePruning_) && (!useKtPruning_) ) {
443  fastjet::Pruner * pruner = new fastjet::Pruner(fastjet::cambridge_algorithm, zCut_, RcutFactor_);
444  transformers.push_back( transformer_ptr(pruner ));
445  }
446 
447  if ( useDynamicFiltering_ ){
448  fastjet::Filter * filter = new fastjet::Filter( fastjet::Filter(&*rFiltDynamic_, fastjet::SelectorNHardest(nFilt_)));
449  transformers.push_back( transformer_ptr(filter));
450  }
451 
452  if ( useKtPruning_ ) {
453  fastjet::Pruner * pruner = new fastjet::Pruner(fastjet::kt_algorithm, zCut_, RcutFactor_);
454  transformers.push_back( transformer_ptr(pruner ));
455  }
456 
457  if ( useSoftDrop_ ) {
458  fastjet::contrib::SoftDrop * sd = new fastjet::contrib::SoftDrop(beta_, zCut_ );
459  transformers.push_back( transformer_ptr(sd) );
460  }
461 
462 
463  for ( std::vector<fastjet::PseudoJet>::const_iterator ijet = tempJets.begin(),
464  ijetEnd = tempJets.end(); ijet != ijetEnd; ++ijet ) {
465 
466  fastjet::PseudoJet transformedJet = *ijet;
467  bool passed = true;
468  for ( transformer_coll::const_iterator itransf = transformers.begin(),
469  itransfEnd = transformers.end(); itransf != itransfEnd; ++itransf ) {
470  if ( transformedJet != 0 ) {
471  transformedJet = (**itransf)(transformedJet);
472  } else {
473  passed=false;
474  }
475  }
476 
477  if ( passed ) {
478  fjJets_.push_back( transformedJet );
479  }
480  }
481  }
482 
483 }
484 
485 
486 
488 // define as cmssw plugin
490 
492 
#define LogDebug(id)
double muMin_
for CMSBoostedTauSeedingAlgorithm : subjet pt min
JetType::Type jetTypeE
T getParameter(std::string const &) const
double yMin_
for CMSBoostedTauSeedingAlgorithm : max mass-drop
bool useFiltering_
Mass-drop tagging for boosted Higgs.
int i
Definition: DBlmapReader.cc:9
double rFilt_
for mass-drop tagging, symmetry cut: min(pt1^2,pt2^2) * dR(1,2) / mjet &gt; ycut
tuple Filter
Definition: Filter_cff.py:5
reco::Particle::Point vertex_
double dRMax_
for CMSBoostedTauSeedingAlgorithm : min dR
double subjetPtMin_
for soft drop : beta (angular exponent)
virtual std::vector< reco::CandidatePtr > getConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents)
std::vector< fastjet::PseudoJet > fjJets_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
virtual void setPileup(float fEnergy)
Set pileup energy contribution as calculated by algorithm.
Definition: Jet.h:108
virtual void inputTowers()
virtual void runAlgorithm(edm::Event &iEvent, const edm::EventSetup &iSetup)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool useDynamicFiltering_
Jet filtering technique.
bool exists(std::string const &parameterName) const
checks if a parameter exists
double RcutFactor_
for pruning OR soft drop: constituent minimum pt fraction of parent cluster
bool useSoftDrop_
constituent subtraction technique
double yMax_
for CMSBoostedTauSeedingAlgorithm : min asymmetry
bool useKtPruning_
algorithm for seeding reconstruction of boosted Taus (similar to mass-drop tagging) ...
FastjetJetProducer(const edm::ParameterSet &iConfig)
bool useTrimming_
Use dynamic filtering radius (as in arXiv:0802.2470)
double csRParam_
for constituent subtraction : maximum rapidity for ghosts
bool makeTrackJet(const JetType::Type &fTag)
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:103
std::vector< fastjet::PseudoJet > fjInputs_
std::vector< transformer_ptr > transformer_coll
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
double trimPtFracMin_
for filtering, pruning: number of subjets expected
int iEvent
Definition: GenABIO.cc:230
bool usePruning_
Jet trimming technique.
std::vector< edm::Ptr< reco::Candidate > > inputs_
bool useCMSBoostedTauSeedingAlgorithm_
Jet pruning technique.
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
ClusterSequencePtr fjClusterSeq_
virtual void produceTrackJets(edm::Event &iEvent, const edm::EventSetup &iSetup)
vector< PseudoJet > jets
edm::EDGetTokenT< edm::View< reco::RecoChargedRefCandidate > > input_chrefcand_token_
for CMSBoostedTauSeedingAlgorithm : max depth for descending into clustering sequence ...
double csRho_EtaMax_
for pruning: constituent dR * pt/2m &lt; rcut_factor
double muMax_
for CMSBoostedTauSeedingAlgorithm : min mass-drop
double muCut_
Soft drop.
double rFiltFactor_
for dynamic filtering radius (as in arXiv:0802.2470)
DynamicRfiltPtr rFiltDynamic_
for filtering, trimming: dR scale of sub-clustering
void setPrimaryVertex(const reco::VertexRef &vtx)
set associated primary vertex
Definition: TrackJet.cc:80
double dRMin_
for CMSBoostedTauSeedingAlgorithm : max asymmetry
boost::shared_ptr< DynamicRfilt > DynamicRfiltPtr
double beta_
for constituent subtraction : R parameter for KT alg in jet median background estimator ...
double zCut_
for trimming: constituent minimum pt fraction of full jet
virtual void setVertex(const Point &vertex)
set vertex
Jets made out of tracks.
Definition: TrackJet.h:27
int nFilt_
for dynamic filtering radius (as in arXiv:0802.2470)
int maxDepth_
for CMSBoostedTauSeedingAlgorithm : max dR
double sd
edm::EDGetTokenT< reco::VertexCollection > input_vertex_token_
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
std::unique_ptr< transformer > transformer_ptr
double yCut_
for mass-drop tagging, m0/mjet (m0 = mass of highest mass subjet)
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
bool useConstituentSubtraction_
Use Kt clustering algorithm for pruning (default is Cambridge/Aachen)
AreaDefinitionPtr fjAreaDefinition_
volatile std::atomic< bool > shutdown_flag false
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, edm::EventSetup const &c)
Definition: JetSpecific.cc:41