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  R0_(-1.0)
80 {
81 
82  if ( iConfig.exists("UseOnlyVertexTracks") )
83  useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
84  else
85  useOnlyVertexTracks_ = false;
86 
87  if ( iConfig.exists("UseOnlyOnePV") )
88  useOnlyOnePV_ = iConfig.getParameter<bool>("UseOnlyOnePV");
89  else
90  useOnlyOnePV_ = false;
91 
92  if ( iConfig.exists("DzTrVtxMax") )
93  dzTrVtxMax_ = iConfig.getParameter<double>("DzTrVtxMax");
94  else
95  dzTrVtxMax_ = 999999.;
96  if ( iConfig.exists("DxyTrVtxMax") )
97  dxyTrVtxMax_ = iConfig.getParameter<double>("DxyTrVtxMax");
98  else
99  dxyTrVtxMax_ = 999999.;
100  if ( iConfig.exists("MinVtxNdof") )
101  minVtxNdof_ = iConfig.getParameter<int>("MinVtxNdof");
102  else
103  minVtxNdof_ = 5;
104  if ( iConfig.exists("MaxVtxZ") )
105  maxVtxZ_ = iConfig.getParameter<double>("MaxVtxZ");
106  else
107  maxVtxZ_ = 15;
108 
109 
110  if ( iConfig.exists("useFiltering") ||
111  iConfig.exists("useTrimming") ||
112  iConfig.exists("usePruning") ||
113  iConfig.exists("useMassDropTagger") ||
114  iConfig.exists("useCMSBoostedTauSeedingAlgorithm") ||
115  iConfig.exists("useConstituentSubtraction") ||
116  iConfig.exists("useSoftDrop")
117  ) {
118  useMassDropTagger_=false;
119  useFiltering_=false;
120  useDynamicFiltering_=false;
121  useTrimming_=false;
122  usePruning_=false;
124  useKtPruning_=false;
126  useSoftDrop_ = false;
127  rFilt_=-1.0;
128  rFiltFactor_=-1.0;
129  nFilt_=-1;
130  trimPtFracMin_=-1.0;
131  zCut_=-1.0;
132  RcutFactor_=-1.0;
133  muCut_=-1.0;
134  yCut_=-1.0;
135  subjetPtMin_ = -1.0;
136  muMin_ = -1.0;
137  muMax_ = -1.0;
138  yMin_ = -1.0;
139  yMax_ = -1.0;
140  dRMin_ = -1.0;
141  dRMax_ = -1.0;
142  maxDepth_ = -1;
143  csRho_EtaMax_ = -1.0;
144  csRParam_ = -1.0;
145  beta_ = -1.0;
146  R0_ = -1.0;
147  useExplicitGhosts_ = true;
148 
149  if ( iConfig.exists("useMassDropTagger") ) {
150  useMassDropTagger_ = iConfig.getParameter<bool>("useMassDropTagger");
151  muCut_ = iConfig.getParameter<double>("muCut");
152  yCut_ = iConfig.getParameter<double>("yCut");
153  }
154 
155  if ( iConfig.exists("useFiltering") ) {
156  useFiltering_ = iConfig.getParameter<bool>("useFiltering");
157  rFilt_ = iConfig.getParameter<double>("rFilt");
158  nFilt_ = iConfig.getParameter<int>("nFilt");
159  if ( iConfig.exists("useDynamicFiltering") ) {
160  useDynamicFiltering_ = iConfig.getParameter<bool>("useDynamicFiltering");
161  rFiltFactor_ = iConfig.getParameter<double>("rFiltFactor");
162  if ( useDynamicFiltering_ )
163  rFiltDynamic_ = DynamicRfiltPtr(new DynamicRfilt(rFilt_, rFiltFactor_));
164  }
165  }
166 
167  if ( iConfig.exists("useTrimming") ) {
168  useTrimming_ = iConfig.getParameter<bool>("useTrimming");
169  rFilt_ = iConfig.getParameter<double>("rFilt");
170  trimPtFracMin_ = iConfig.getParameter<double>("trimPtFracMin");
171  }
172 
173  if ( iConfig.exists("usePruning") ) {
174  usePruning_ = iConfig.getParameter<bool>("usePruning");
175  zCut_ = iConfig.getParameter<double>("zcut");
176  RcutFactor_ = iConfig.getParameter<double>("rcut_factor");
177  nFilt_ = iConfig.getParameter<int>("nFilt");
178  if ( iConfig.exists("useKtPruning") )
179  useKtPruning_ = iConfig.getParameter<bool>("useKtPruning");
180  }
181 
182  if ( iConfig.exists("useCMSBoostedTauSeedingAlgorithm") ) {
183  useCMSBoostedTauSeedingAlgorithm_ = iConfig.getParameter<bool>("useCMSBoostedTauSeedingAlgorithm");
184  subjetPtMin_ = iConfig.getParameter<double>("subjetPtMin");
185  muMin_ = iConfig.getParameter<double>("muMin");
186  muMax_ = iConfig.getParameter<double>("muMax");
187  yMin_ = iConfig.getParameter<double>("yMin");
188  yMax_ = iConfig.getParameter<double>("yMax");
189  dRMin_ = iConfig.getParameter<double>("dRMin");
190  dRMax_ = iConfig.getParameter<double>("dRMax");
191  maxDepth_ = iConfig.getParameter<int>("maxDepth");
192  }
193 
194  if ( iConfig.exists("useConstituentSubtraction") ) {
195 
196  if ( fjAreaDefinition_.get() == 0 ) {
197  throw cms::Exception("AreaMustBeSet") << "Logic error. The area definition must be set if you use constituent subtraction." << std::endl;
198  }
199 
200  useConstituentSubtraction_ = iConfig.getParameter<bool>("useConstituentSubtraction");
201  csRho_EtaMax_ = iConfig.getParameter<double>("csRho_EtaMax");
202  csRParam_ = iConfig.getParameter<double>("csRParam");
203  }
204 
205  if ( iConfig.exists("useSoftDrop") ) {
206  if ( usePruning_ ) {
207  throw cms::Exception("PruningAndSoftDrop") << "Logic error. Soft drop is a generalized pruning, do not run them together." << std::endl;
208  }
209  useSoftDrop_ = iConfig.getParameter<bool>("useSoftDrop");
210  zCut_ = iConfig.getParameter<double>("zcut");
211  beta_ = iConfig.getParameter<double>("beta");
212  R0_ = iConfig.getParameter<double>("R0");
213  }
214 
215  }
216 
217  input_chrefcand_token_ = consumes<edm::View<reco::RecoChargedRefCandidate> >(src_);
218 
219 }
220 
221 
222 //______________________________________________________________________________
224 {
225 }
226 
227 
229 // implementation of member functions
231 
233 {
234 
235  // for everything but track jets
236  if (!makeTrackJet(jetTypeE)) {
237 
238  // use the default production from one collection
239  VirtualJetProducer::produce( iEvent, iSetup );
240 
241  } else { // produce trackjets from tracks grouped per primary vertex
242 
243  produceTrackJets(iEvent, iSetup);
244 
245  }
246 
247 }
248 
249 
251 {
252 
253  // read in the track candidates
255  iEvent.getByToken(input_chrefcand_token_, inputsHandle);
256 
257  // make collection with pointers so we can play around with it
258  std::vector<edm::Ptr<reco::RecoChargedRefCandidate> > allInputs;
259  std::vector<edm::Ptr<reco::Candidate> > origInputs;
260  for (size_t i = 0; i < inputsHandle->size(); ++i) {
261  allInputs.push_back(inputsHandle->ptrAt(i));
262  origInputs.push_back(inputsHandle->ptrAt(i));
263  }
264 
265  // read in the PV collection
267  iEvent.getByToken(input_vertex_token_, pvCollection);
268  // define the overall output jet container
269  std::auto_ptr<std::vector<reco::TrackJet> > jets(new std::vector<reco::TrackJet>() );
270 
271  // loop over the good vertices, clustering for each vertex separately
272  for (reco::VertexCollection::const_iterator itVtx = pvCollection->begin(); itVtx != pvCollection->end(); ++itVtx) {
273  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
274 
275  // clear the intermediate containers
276  inputs_.clear();
277  fjInputs_.clear();
278  fjJets_.clear();
279 
280  // if only vertex-associated tracks should be used
281  if (useOnlyVertexTracks_) {
282  // loop over the tracks associated to the vertex
283  for (reco::Vertex::trackRef_iterator itTr = itVtx->tracks_begin(); itTr != itVtx->tracks_end(); ++itTr) {
284  // whether a match was found in the track candidate input
285  bool found = false;
286  // loop over input track candidates
287  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
288  // match the input track candidate to the track from the vertex
289  reco::TrackRef trref(itTr->castTo<reco::TrackRef>());
290  // check if the tracks match
291  if ((*itIn)->track() == trref) {
292  found = true;
293  // add this track candidate to the input for clustering
294  inputs_.push_back(*itIn);
295  // erase the track candidate from the total list of input, so we don't reuse it later
296  allInputs.erase(itIn);
297  // found the candidate track corresponding to the vertex track, so stop the loop
298  break;
299  } // end if match found
300  } // end loop over input tracks
301  // give an info message in case no match is found (can happen if candidates are subset of tracks used for clustering)
302  if (!found) edm::LogInfo("FastjetTrackJetProducer") << "Ignoring a track at vertex which is not in input track collection!";
303  } // end loop over tracks associated to vertex
304  // if all inpt track candidates should be used
305  } else {
306  // loop over input track candidates
307  for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
308  // check if the track is close enough to the vertex
309  float dz = (*itIn)->track()->dz(itVtx->position());
310  float dxy = (*itIn)->track()->dxy(itVtx->position());
311  if (fabs(dz) > dzTrVtxMax_) continue;
312  if (fabs(dxy) > dxyTrVtxMax_) continue;
313  bool closervtx = false;
314  // now loop over the good vertices a second time
315  for (reco::VertexCollection::const_iterator itVtx2 = pvCollection->begin(); itVtx2 != pvCollection->end(); ++itVtx2) {
316  if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
317  // and check this track is closer to any other vertex (if more than 1 vertex considered)
318  if (!useOnlyOnePV_ &&
319  itVtx != itVtx2 &&
320  fabs((*itIn)->track()->dz(itVtx2->position())) < fabs(dz)) {
321  closervtx = true;
322  break; // 1 closer vertex makes the track already not matched, so break
323  }
324  }
325  // don't add this track if another vertex is found closer
326  if (closervtx) continue;
327  // add this track candidate to the input for clustering
328  inputs_.push_back(*itIn);
329  // erase the track candidate from the total list of input, so we don't reuse it later
330  allInputs.erase(itIn);
331  // take a step back in the loop since we just erased
332  --itIn;
333  }
334  }
335 
336  // convert candidates in inputs_ to fastjet::PseudoJets in fjInputs_
337  fjInputs_.reserve(inputs_.size());
338  inputTowers();
339  LogDebug("FastjetTrackJetProducer") << "Inputted towers\n";
340 
341  // run algorithm, using fjInputs_, modifying fjJets_ and allocating fjClusterSeq_
342  runAlgorithm(iEvent, iSetup);
343  LogDebug("FastjetTrackJetProducer") << "Ran algorithm\n";
344 
345  // convert our jets and add to the overall jet vector
346  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
347  // get the constituents from fastjet
348  std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJets_[ijet]));
349  // convert them to CandidatePtr vector
350  std::vector<reco::CandidatePtr> constituents = getConstituents(fjConstituents);
351  // fill the trackjet
353  // write the specifics to the jet (simultaneously sets 4-vector, vertex).
354  writeSpecific( jet,
355  reco::Particle::LorentzVector(fjJets_[ijet].px(), fjJets_[ijet].py(), fjJets_[ijet].pz(), fjJets_[ijet].E()),
356  vertex_, constituents, iSetup);
357  jet.setJetArea(0);
358  jet.setPileup(0);
359  jet.setPrimaryVertex(edm::Ref<reco::VertexCollection>(pvCollection, (int) (itVtx-pvCollection->begin())));
360  jet.setVertex(itVtx->position());
361  jets->push_back(jet);
362  }
363 
364  if (useOnlyOnePV_) break; // stop vertex loop if only one vertex asked for
365  } // end loop over vertices
366 
367  // put the jets in the collection
368  LogDebug("FastjetTrackJetProducer") << "Put " << jets->size() << " jets in the event.\n";
369  iEvent.put(jets);
370 
371 }
372 
373 
374 //______________________________________________________________________________
376 {
377  // run algorithm
378  /*
379  fjInputs_.clear();
380  double px, py , pz, E;
381  string line;
382  std::ifstream fin("dump3.txt");
383  while (getline(fin, line)){
384  if (line == "#END") break;
385  if (line.substr(0,1) == "#") {continue;}
386  istringstream istr(line);
387  istr >> px >> py >> pz >> E;
388  // create a fastjet::PseudoJet with these components and put it onto
389  // back of the input_particles vector
390  fastjet::PseudoJet j(px,py,pz,E);
391  //if ( fabs(j.rap()) < inputEtaMax )
392  fjInputs_.push_back(fastjet::PseudoJet(px,py,pz,E));
393  }
394  fin.close();
395  */
396 
397  if ( !doAreaFastjet_ && !doRhoFastjet_) {
398  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
399  } else if (voronoiRfact_ <= 0) {
400  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjAreaDefinition_ ) );
401  } else {
402  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
403  }
404 
406  fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
407  } else {
408  fjJets_.clear();
409 
410 
411  transformer_coll transformers;
412 
413 
414  std::vector<fastjet::PseudoJet> tempJets = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
415 
416  unique_ptr<fastjet::JetMedianBackgroundEstimator> bge_rho;
418  fastjet::Selector rho_range = fastjet::SelectorAbsRapMax(csRho_EtaMax_);
419  bge_rho = unique_ptr<fastjet::JetMedianBackgroundEstimator> (new fastjet::JetMedianBackgroundEstimator(rho_range, fastjet::JetDefinition(fastjet::kt_algorithm, csRParam_), *fjAreaDefinition_) );
420  bge_rho->set_particles(fjInputs_);
421  fastjet::contrib::ConstituentSubtractor * constituentSubtractor = new fastjet::contrib::ConstituentSubtractor(bge_rho.get());
422  // this sets the same background estimator to be used for deltaMass density, rho_m, as for pt density, rho:
423  constituentSubtractor->use_common_bge_for_rho_and_rhom(true); // for massless input particles it does not make any difference (rho_m is always zero)
424 
425  transformers.push_back( transformer_ptr(constituentSubtractor) );
426  };
427  if ( useMassDropTagger_ ) {
428  fastjet::MassDropTagger * md_tagger = new fastjet::MassDropTagger ( muCut_, yCut_ );
429  transformers.push_back( transformer_ptr(md_tagger) );
430  }
432  fastjet::contrib::CMSBoostedTauSeedingAlgorithm * tau_tagger =
433  new fastjet::contrib::CMSBoostedTauSeedingAlgorithm ( subjetPtMin_, muMin_, muMax_, yMin_, yMax_, dRMin_, dRMax_, maxDepth_, verbosity_ );
434  transformers.push_back( transformer_ptr(tau_tagger ));
435  }
436  if ( useTrimming_ ) {
437  fastjet::Filter * trimmer = new fastjet::Filter(fastjet::JetDefinition(fastjet::kt_algorithm, rFilt_), fastjet::SelectorPtFractionMin(trimPtFracMin_));
438  transformers.push_back( transformer_ptr(trimmer) );
439  }
440  if ( (useFiltering_) && (!useDynamicFiltering_) ) {
441  fastjet::Filter * filter = new fastjet::Filter(fastjet::JetDefinition(fastjet::cambridge_algorithm, rFilt_), fastjet::SelectorNHardest(nFilt_));
442  transformers.push_back( transformer_ptr(filter));
443  }
444 
445  if ( (usePruning_) && (!useKtPruning_) ) {
446  fastjet::Pruner * pruner = new fastjet::Pruner(fastjet::cambridge_algorithm, zCut_, RcutFactor_);
447  transformers.push_back( transformer_ptr(pruner ));
448  }
449 
450  if ( useDynamicFiltering_ ){
451  fastjet::Filter * filter = new fastjet::Filter( fastjet::Filter(&*rFiltDynamic_, fastjet::SelectorNHardest(nFilt_)));
452  transformers.push_back( transformer_ptr(filter));
453  }
454 
455  if ( useKtPruning_ ) {
456  fastjet::Pruner * pruner = new fastjet::Pruner(fastjet::kt_algorithm, zCut_, RcutFactor_);
457  transformers.push_back( transformer_ptr(pruner ));
458  }
459 
460  if ( useSoftDrop_ ) {
461  fastjet::contrib::SoftDrop * sd = new fastjet::contrib::SoftDrop(beta_, zCut_, R0_ );
462  transformers.push_back( transformer_ptr(sd) );
463  }
464 
465 
466  for ( std::vector<fastjet::PseudoJet>::const_iterator ijet = tempJets.begin(),
467  ijetEnd = tempJets.end(); ijet != ijetEnd; ++ijet ) {
468 
469  fastjet::PseudoJet transformedJet = *ijet;
470  bool passed = true;
471  for ( transformer_coll::const_iterator itransf = transformers.begin(),
472  itransfEnd = transformers.end(); itransf != itransfEnd; ++itransf ) {
473  if ( transformedJet != 0 ) {
474  transformedJet = (**itransf)(transformedJet);
475  } else {
476  passed=false;
477  }
478  }
479 
480  if ( passed ) {
481  fjJets_.push_back( transformedJet );
482  }
483  }
484  }
485 
486 }
487 
488 
489 
491 // define as cmssw plugin
493 
495 
#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 : R0 (angular distance normalization - should be set to jet radius in most cases) ...
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:449
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_
double R0_
for soft drop : beta (angular exponent)
std::unique_ptr< transformer > transformer_ptr
double yCut_
for mass-drop tagging, m0/mjet (m0 = mass of highest mass subjet)
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
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:21
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