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