CMS 3D CMS Logo

VirtualJetProducer.cc
Go to the documentation of this file.
1 //
3 // VirtualJetProducer
4 // ------------------
5 //
6 // 04/21/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
8 
14 
18 
27 
42 
45 
46 #include "fastjet/SISConePlugin.hh"
47 #include "fastjet/CMSIterativeConePlugin.hh"
48 #include "fastjet/ATLASConePlugin.hh"
49 #include "fastjet/CDFMidPointPlugin.hh"
50 
51 #include <iostream>
52 #include <memory>
53 #include <algorithm>
54 #include <limits>
55 #include <cmath>
56 #include <vdt/vdtMath.h>
57 
58 using namespace std;
59 using namespace edm;
60 
61 namespace reco {
62  namespace helper {
64  bool operator()(const fastjet::PseudoJet& t1, const fastjet::PseudoJet& t2) const {
65  return t1.perp2() > t2.perp2();
66  }
67  };
68 
69  } // namespace helper
70 } // namespace reco
71 
72 //______________________________________________________________________________
73 const char* const VirtualJetProducer::JetType::names[] = {
74  "BasicJet", "GenJet", "CaloJet", "PFJet", "TrackJet", "PFClusterJet"};
75 
76 //______________________________________________________________________________
78  const char* const* pos = std::find(names, names + LastJetType, name);
79  if (pos == names + LastJetType) {
80  std::string errorMessage = "Requested jetType not supported: " + name + "\n";
81  throw cms::Exception("Configuration", errorMessage);
82  }
83  return (Type)(pos - names);
84 }
85 
87  if (writeCompound_) {
88  produces<reco::BasicJetCollection>();
89  }
90 
91  if (writeJetsWithConst_) {
92  produces<reco::PFCandidateCollection>(tag).setBranchAlias(alias);
93  produces<reco::PFJetCollection>();
94  } else {
95  if (makeCaloJet(jetTypeE)) {
96  produces<reco::CaloJetCollection>(tag).setBranchAlias(alias);
97  } else if (makePFJet(jetTypeE)) {
98  produces<reco::PFJetCollection>(tag).setBranchAlias(alias);
99  } else if (makeGenJet(jetTypeE)) {
100  produces<reco::GenJetCollection>(tag).setBranchAlias(alias);
101  } else if (makeTrackJet(jetTypeE)) {
102  produces<reco::TrackJetCollection>(tag).setBranchAlias(alias);
103  } else if (makePFClusterJet(jetTypeE)) {
104  produces<reco::PFClusterJetCollection>(tag).setBranchAlias(alias);
105  } else if (makeBasicJet(jetTypeE)) {
106  produces<reco::BasicJetCollection>(tag).setBranchAlias(alias);
107  }
108  }
109 }
110 
112 // construction / destruction
114 
115 //______________________________________________________________________________
117  moduleLabel_ = iConfig.getParameter<string>("@module_label");
118  src_ = iConfig.getParameter<edm::InputTag>("src");
119  srcPVs_ = iConfig.getParameter<edm::InputTag>("srcPVs");
120  jetType_ = iConfig.getParameter<string>("jetType");
121  jetAlgorithm_ = iConfig.getParameter<string>("jetAlgorithm");
122  rParam_ = iConfig.getParameter<double>("rParam");
123  inputEtMin_ = iConfig.getParameter<double>("inputEtMin");
124  inputEMin_ = iConfig.getParameter<double>("inputEMin");
125  jetPtMin_ = iConfig.getParameter<double>("jetPtMin");
126  doPVCorrection_ = iConfig.getParameter<bool>("doPVCorrection");
127  doAreaFastjet_ = iConfig.getParameter<bool>("doAreaFastjet");
128  doRhoFastjet_ = iConfig.getParameter<bool>("doRhoFastjet");
129  jetCollInstanceName_ = iConfig.getParameter<string>("jetCollInstanceName");
130  doPUOffsetCorr_ = iConfig.getParameter<bool>("doPUOffsetCorr");
131  puSubtractorName_ = iConfig.getParameter<string>("subtractorName");
132  useExplicitGhosts_ =
133  iConfig.getParameter<bool>("useExplicitGhosts"); // use explicit ghosts in the fastjet clustering sequence?
134  doAreaDiskApprox_ = iConfig.getParameter<bool>("doAreaDiskApprox");
135  voronoiRfact_ = iConfig.getParameter<double>(
136  "voronoiRfact"); // Voronoi-based area calculation allows for an empirical scale factor
137  rhoEtaMax_ = iConfig.getParameter<double>(
138  "Rho_EtaMax"); // do fasjet area / rho calcluation? => accept corresponding parameters
139  ghostEtaMax_ = iConfig.getParameter<double>("Ghost_EtaMax");
140  activeAreaRepeats_ = iConfig.getParameter<int>("Active_Area_Repeats");
141  ghostArea_ = iConfig.getParameter<double>("GhostArea");
142  restrictInputs_ = iConfig.getParameter<bool>("restrictInputs"); // restrict inputs to first "maxInputs" towers?
143  maxInputs_ = iConfig.getParameter<unsigned int>("maxInputs");
144  writeCompound_ = iConfig.getParameter<bool>(
145  "writeCompound"); // Check to see if we are writing compound jets for substructure and jet grooming
146  writeJetsWithConst_ = iConfig.getParameter<bool>("writeJetsWithConst"); //write subtracted jet constituents
147  doFastJetNonUniform_ = iConfig.getParameter<bool>("doFastJetNonUniform");
148  puCenters_ = iConfig.getParameter<vector<double>>("puCenters");
149  puWidth_ = iConfig.getParameter<double>("puWidth");
150  nExclude_ = iConfig.getParameter<unsigned int>("nExclude");
151  useDeterministicSeed_ = iConfig.getParameter<bool>("useDeterministicSeed");
152  minSeed_ = iConfig.getParameter<unsigned int>("minSeed");
153  verbosity_ = iConfig.getParameter<int>("verbosity");
154  applyWeight_ = iConfig.getParameter<bool>("applyWeight");
155  if (applyWeight_) {
156  edm::InputTag srcWeights = iConfig.getParameter<edm::InputTag>("srcWeights");
157  if (srcWeights.label() == src_.label())
158  LogWarning("VirtualJetProducer")
159  << "Particle and weights collection have the same label. You may be applying the same weights twice.\n";
160  if (!srcWeights.label().empty())
161  input_weights_token_ = consumes<edm::ValueMap<float>>(srcWeights);
162  }
163 
164  anomalousTowerDef_ = std::make_unique<AnomalousTower>(iConfig);
165 
166  input_vertex_token_ = consumes<reco::VertexCollection>(srcPVs_);
167  input_candidateview_token_ = consumes<reco::CandidateView>(src_);
168  input_candidatefwdptr_token_ =
169  consumes<vector<edm::FwdPtr<reco::PFCandidate>>>(iConfig.getParameter<edm::InputTag>("src"));
170  input_packedcandidatefwdptr_token_ =
171  consumes<vector<edm::FwdPtr<pat::PackedCandidate>>>(iConfig.getParameter<edm::InputTag>("src"));
172  input_gencandidatefwdptr_token_ =
173  consumes<vector<edm::FwdPtr<reco::GenParticle>>>(iConfig.getParameter<edm::InputTag>("src"));
174  input_packedgencandidatefwdptr_token_ =
175  consumes<vector<edm::FwdPtr<pat::PackedGenParticle>>>(iConfig.getParameter<edm::InputTag>("src"));
176 
177  //
178  // additional parameters to think about:
179  // - overlap threshold (set to 0.75 for the time being)
180  // - p parameter for generalized kT (set to -2 for the time being)
181  // - fastjet PU subtraction parameters (not yet considered)
182  //
183  if (jetAlgorithm_ == "Kt")
184  fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(fastjet::kt_algorithm, rParam_);
185 
186  else if (jetAlgorithm_ == "CambridgeAachen")
187  fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(fastjet::cambridge_algorithm, rParam_);
188 
189  else if (jetAlgorithm_ == "AntiKt")
190  fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(fastjet::antikt_algorithm, rParam_);
191 
192  else if (jetAlgorithm_ == "GeneralizedKt")
193  fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(fastjet::genkt_algorithm, rParam_, -2);
194 
195  else if (jetAlgorithm_ == "SISCone") {
196  fjPlugin_ = PluginPtr(new fastjet::SISConePlugin(rParam_, 0.75, 0, 0.0, false, fastjet::SISConePlugin::SM_pttilde));
197  fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(&*fjPlugin_);
198 
199  } else if (jetAlgorithm_ == "IterativeCone") {
200  fjPlugin_ = PluginPtr(new fastjet::CMSIterativeConePlugin(rParam_, 1.0));
201  fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(&*fjPlugin_);
202 
203  } else if (jetAlgorithm_ == "CDFMidPoint") {
204  fjPlugin_ = PluginPtr(new fastjet::CDFMidPointPlugin(rParam_, 0.75));
205  fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(&*fjPlugin_);
206 
207  } else if (jetAlgorithm_ == "ATLASCone") {
208  fjPlugin_ = PluginPtr(new fastjet::ATLASConePlugin(rParam_));
209  fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(&*fjPlugin_);
210 
211  } else {
212  throw cms::Exception("Invalid jetAlgorithm") << "Jet algorithm for VirtualJetProducer is invalid, Abort!\n";
213  }
214 
215  jetTypeE = JetType::byName(jetType_);
216 
217  if (JetType::CaloJet == jetTypeE) {
218  geometry_token_ = esConsumes();
219  topology_token_ = esConsumes();
220  }
221 
222  if (doPUOffsetCorr_) {
223  if (puSubtractorName_.empty()) {
224  LogWarning("VirtualJetProducer")
225  << "Pile Up correction on; however, pile up type is not specified. Using default... \n";
226  subtractor_ = std::make_shared<PileUpSubtractor>(iConfig, consumesCollector());
227  } else
228  subtractor_ = std::shared_ptr<PileUpSubtractor>(
229  PileUpSubtractorFactory::get()->create(puSubtractorName_, iConfig, consumesCollector()));
230  }
231 
232  // do approximate disk-based area calculation => warn if conflicting request
233  if (doAreaDiskApprox_ && doAreaFastjet_)
234  throw cms::Exception("Conflicting area calculations")
235  << "Both the calculation of jet area via fastjet and via an analytical disk approximation have been requested. "
236  "Please decide on one.\n";
237 
238  if (doAreaFastjet_ || doRhoFastjet_) {
239  if (voronoiRfact_ <= 0) {
240  fjActiveArea_ = std::make_shared<fastjet::GhostedAreaSpec>(ghostEtaMax_, activeAreaRepeats_, ghostArea_);
241 
242  if (!useExplicitGhosts_) {
243  fjAreaDefinition_ = std::make_shared<fastjet::AreaDefinition>(fastjet::active_area, *fjActiveArea_);
244  } else {
245  fjAreaDefinition_ =
246  std::make_shared<fastjet::AreaDefinition>(fastjet::active_area_explicit_ghosts, *fjActiveArea_);
247  }
248  }
249  fjSelector_ = std::make_shared<fastjet::Selector>(fastjet::SelectorAbsRapMax(rhoEtaMax_));
250  }
251 
252  if ((doFastJetNonUniform_) && (puCenters_.empty()))
253  throw cms::Exception("doFastJetNonUniform")
254  << "Parameter puCenters for doFastJetNonUniform is not defined." << std::endl;
255 
256  // make the "produces" statements
257  makeProduces(moduleLabel_, jetCollInstanceName_);
258  produces<vector<double>>("rhos");
259  produces<vector<double>>("sigmas");
260  produces<double>("rho");
261  produces<double>("sigma");
262 }
263 
264 //______________________________________________________________________________
266 
268 // implementation of member functions
270 
271 //______________________________________________________________________________
273  // If requested, set the fastjet random seed to a deterministic function
274  // of the run/lumi/event.
275  // NOTE!!! The fastjet random number sequence is a global singleton.
276  // Thus, we have to create an object and get access to the global singleton
277  // in order to change it.
278  if (useDeterministicSeed_) {
279  fastjet::GhostedAreaSpec gas;
280  std::vector<int> seeds(2);
281  unsigned int runNum_uint = static_cast<unsigned int>(iEvent.id().run());
282  unsigned int evNum_uint = static_cast<unsigned int>(iEvent.id().event());
283  seeds[0] = std::max(runNum_uint, minSeed_ + 3) + 3 * evNum_uint;
284  seeds[1] = std::max(runNum_uint, minSeed_ + 5) + 5 * evNum_uint;
285  gas.set_random_status(seeds);
286  }
287 
288  LogDebug("VirtualJetProducer") << "Entered produce\n";
289  //determine signal vertex2
290  vertex_ = reco::Jet::Point(0, 0, 0);
291  if ((makeCaloJet(jetTypeE) || makePFJet(jetTypeE)) && doPVCorrection_) {
292  LogDebug("VirtualJetProducer") << "Adding PV info\n";
294  iEvent.getByToken(input_vertex_token_, pvCollection);
295  if (!pvCollection->empty())
296  vertex_ = pvCollection->begin()->position();
297  }
298 
299  // Get Weights Collection
300  if ((applyWeight_) && (!input_weights_token_.isUninitialized()))
301  weights_ = iEvent.get(input_weights_token_);
302 
303  // For Pileup subtraction using offset correction:
304  // set up geometry map
305  if (doPUOffsetCorr_) {
306  subtractor_->setupGeometryMap(iEvent, iSetup);
307  }
308 
309  // clear data
310  LogDebug("VirtualJetProducer") << "Clear data\n";
311  fjInputs_.clear();
312  fjJets_.clear();
313  inputs_.clear();
314 
315  // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
317 
321  edm::Handle<std::vector<edm::FwdPtr<pat::PackedGenParticle>>> packedgeninputsHandleAsFwdPtr;
322 
323  bool isView = iEvent.getByToken(input_candidateview_token_, inputsHandle);
324  if (isView) {
325  if (inputsHandle->empty()) {
326  output(iEvent, iSetup);
327  return;
328  }
329  for (size_t i = 0; i < inputsHandle->size(); ++i) {
330  inputs_.push_back(inputsHandle->ptrAt(i));
331  }
332  } else {
333  bool isPF = iEvent.getByToken(input_candidatefwdptr_token_, pfinputsHandleAsFwdPtr);
334  bool isPFFwdPtr = iEvent.getByToken(input_packedcandidatefwdptr_token_, packedinputsHandleAsFwdPtr);
335  bool isGen = iEvent.getByToken(input_gencandidatefwdptr_token_, geninputsHandleAsFwdPtr);
336  bool isGenFwdPtr = iEvent.getByToken(input_packedgencandidatefwdptr_token_, packedgeninputsHandleAsFwdPtr);
337 
338  if (isPF) {
339  if (pfinputsHandleAsFwdPtr->empty()) {
340  output(iEvent, iSetup);
341  return;
342  }
343  for (size_t i = 0; i < pfinputsHandleAsFwdPtr->size(); ++i) {
344  if ((*pfinputsHandleAsFwdPtr)[i].ptr().isAvailable()) {
345  inputs_.push_back((*pfinputsHandleAsFwdPtr)[i].ptr());
346  } else if ((*pfinputsHandleAsFwdPtr)[i].backPtr().isAvailable()) {
347  inputs_.push_back((*pfinputsHandleAsFwdPtr)[i].backPtr());
348  }
349  }
350  } else if (isPFFwdPtr) {
351  if (packedinputsHandleAsFwdPtr->empty()) {
352  output(iEvent, iSetup);
353  return;
354  }
355  for (size_t i = 0; i < packedinputsHandleAsFwdPtr->size(); ++i) {
356  if ((*packedinputsHandleAsFwdPtr)[i].ptr().isAvailable()) {
357  inputs_.push_back((*packedinputsHandleAsFwdPtr)[i].ptr());
358  } else if ((*packedinputsHandleAsFwdPtr)[i].backPtr().isAvailable()) {
359  inputs_.push_back((*packedinputsHandleAsFwdPtr)[i].backPtr());
360  }
361  }
362  } else if (isGen) {
363  if (geninputsHandleAsFwdPtr->empty()) {
364  output(iEvent, iSetup);
365  return;
366  }
367  for (size_t i = 0; i < geninputsHandleAsFwdPtr->size(); ++i) {
368  if ((*geninputsHandleAsFwdPtr)[i].ptr().isAvailable()) {
369  inputs_.push_back((*geninputsHandleAsFwdPtr)[i].ptr());
370  } else if ((*geninputsHandleAsFwdPtr)[i].backPtr().isAvailable()) {
371  inputs_.push_back((*geninputsHandleAsFwdPtr)[i].backPtr());
372  }
373  }
374  } else if (isGenFwdPtr) {
375  if (geninputsHandleAsFwdPtr->empty()) {
376  output(iEvent, iSetup);
377  return;
378  }
379  for (size_t i = 0; i < packedgeninputsHandleAsFwdPtr->size(); ++i) {
380  if ((*packedgeninputsHandleAsFwdPtr)[i].ptr().isAvailable()) {
381  inputs_.push_back((*packedgeninputsHandleAsFwdPtr)[i].ptr());
382  } else if ((*packedgeninputsHandleAsFwdPtr)[i].backPtr().isAvailable()) {
383  inputs_.push_back((*packedgeninputsHandleAsFwdPtr)[i].backPtr());
384  }
385  }
386  } else {
387  throw cms::Exception("Invalid Jet Inputs")
388  << "Did not specify appropriate inputs for VirtualJetProducer, Abort!\n";
389  }
390  }
391  LogDebug("VirtualJetProducer") << "Got inputs\n";
392 
393  // Convert candidates to fastjet::PseudoJets.
394  // Also correct to Primary Vertex. Will modify fjInputs_
395  // and use inputs_
396  fjInputs_.reserve(inputs_.size());
397  inputTowers();
398  LogDebug("VirtualJetProducer") << "Inputted towers\n";
399 
400  // For Pileup subtraction using offset correction:
401  // Subtract pedestal.
402  if (doPUOffsetCorr_) {
403  subtractor_->setDefinition(fjJetDefinition_);
404  subtractor_->reset(inputs_, fjInputs_, fjJets_);
405  subtractor_->calculatePedestal(fjInputs_);
406  subtractor_->subtractPedestal(fjInputs_);
407  LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
408  }
409  // Run algorithm. Will modify fjJets_ and allocate fjClusterSeq_.
410  // This will use fjInputs_
411  runAlgorithm(iEvent, iSetup);
412 
413  // if ( doPUOffsetCorr_ ) {
414  // subtractor_->setAlgorithm(fjClusterSeq_);
415  // }
416 
417  LogDebug("VirtualJetProducer") << "Ran algorithm\n";
418  // For Pileup subtraction using offset correction:
419  // Now we find jets and need to recalculate their energy,
420  // mark towers participated in jet,
421  // remove occupied towers from the list and recalculate mean and sigma
422  // put the initial towers collection to the jet,
423  // and subtract from initial towers in jet recalculated mean and sigma of towers
424  if (doPUOffsetCorr_) {
425  LogDebug("VirtualJetProducer") << "Do PUOffsetCorr\n";
426  vector<fastjet::PseudoJet> orphanInput;
427  subtractor_->calculateOrphanInput(orphanInput);
428  subtractor_->calculatePedestal(orphanInput);
429  subtractor_->offsetCorrectJets();
430  }
431  // Write the output jets.
432  // This will (by default) call the member function template
433  // "writeJets", but can be overridden.
434  // this will use inputs_
435  output(iEvent, iSetup);
436  LogDebug("VirtualJetProducer") << "Wrote jets\n";
437 
438  // Clear the work vectors so that memory is free for other modules.
439  // Use the trick of swapping with an empty vector so that the memory
440  // is actually given back rather than silently kept.
441  decltype(fjInputs_)().swap(fjInputs_);
442  decltype(fjJets_)().swap(fjJets_);
443  decltype(inputs_)().swap(inputs_);
444 
445  return;
446 }
447 
448 //______________________________________________________________________________
449 
451  auto inBegin = inputs_.begin(), inEnd = inputs_.end(), i = inBegin;
452  for (; i != inEnd; ++i) {
453  auto const& input = **i;
454  // std::cout << "CaloTowerVI jets " << input->pt() << " " << input->et() << ' '<< input->energy() << ' ' << (isAnomalousTower(input) ? " bad" : " ok") << std::endl;
455  if (edm::isNotFinite(input.pt()))
456  continue;
457  if (input.et() < inputEtMin_)
458  continue;
459  if (input.energy() < inputEMin_)
460  continue;
461  if (isAnomalousTower(*i))
462  continue;
463  // Change by SRR : this is no longer an error nor warning, this can happen with PU mitigation algos.
464  // Also switch to something more numerically safe. (VI: 10^-42GeV????)
465  if (input.pt() < 100 * std::numeric_limits<double>::epsilon()) {
466  continue;
467  }
468  if (makeCaloJet(jetTypeE) && doPVCorrection_) {
469  const CaloTower& tower = dynamic_cast<const CaloTower&>(input);
470  auto const& ct = tower.p4(vertex_); // very expensive as computed in eta/phi
471  fjInputs_.emplace_back(ct.px(), ct.py(), ct.pz(), ct.energy());
472  fjInputs_.back().set_user_index(i - inBegin);
473  //std::cout << "tower:" << *tower << '\n';
474  } else {
475  /*
476  if(makePFJet(jetTypeE)) {
477  reco::PFCandidate& pfc = (reco::PFCandidate&)input;
478  std::cout << "PF cand:" << pfc << '\n';
479  }
480  */
481  if (!applyWeight_) {
482  fjInputs_.emplace_back(input.px(), input.py(), input.pz(), input.energy());
483  fjInputs_.back().set_user_index(i - inBegin);
484  } else {
485  if (input_weights_token_.isUninitialized())
486  throw cms::Exception("InvalidInput")
487  << "applyWeight set to True, but no weights given in VirtualJetProducer\n";
488  float w = weights_[*i];
489  if (w > 0) {
490  fjInputs_.emplace_back(input.px() * w, input.py() * w, input.pz() * w, input.energy() * w);
491  fjInputs_.back().set_user_index(i - inBegin);
492  }
493  }
494  }
495  }
496 
497  if (restrictInputs_ && fjInputs_.size() > maxInputs_) {
499  std::sort(fjInputs_.begin(), fjInputs_.end(), pTComparator);
500  fjInputs_.resize(maxInputs_);
501  edm::LogWarning("JetRecoTooManyEntries")
502  << "Too many inputs in the event, limiting to first " << maxInputs_ << ". Output is suspect.";
503  }
504 }
505 
506 //______________________________________________________________________________
508  if (!makeCaloJet(jetTypeE))
509  return false;
510  else
511  return (*anomalousTowerDef_)(*input);
512 }
513 
514 //------------------------------------------------------------------------------
515 // This is pure virtual.
516 //______________________________________________________________________________
517 // void VirtualJetProducer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup,
518 // std::vector<edm::Ptr<reco::Candidate> > const & inputs_);
519 
520 //______________________________________________________________________________
521 void VirtualJetProducer::copyConstituents(const vector<fastjet::PseudoJet>& fjConstituents, reco::Jet* jet) {
522  for (unsigned int i = 0; i < fjConstituents.size(); ++i) {
523  int index = fjConstituents[i].user_index();
524  if (index >= 0 && static_cast<unsigned int>(index) < inputs_.size())
525  jet->addDaughter(inputs_[index]);
526  }
527 }
528 
529 //______________________________________________________________________________
530 vector<reco::CandidatePtr> VirtualJetProducer::getConstituents(const vector<fastjet::PseudoJet>& fjConstituents) {
531  vector<reco::CandidatePtr> result;
532  result.reserve(fjConstituents.size() / 2);
533  for (unsigned int i = 0; i < fjConstituents.size(); i++) {
534  auto index = fjConstituents[i].user_index();
535  if (index >= 0 && static_cast<unsigned int>(index) < inputs_.size()) {
536  result.emplace_back(inputs_[index]);
537  }
538  }
539  return result;
540 }
541 
542 //_____________________________________________________________________________
543 
545  // Write jets and constitutents. Will use fjJets_, inputs_
546  // and fjClusterSeq_
547 
548  if (writeCompound_) {
549  // Write jets and subjets
550  switch (jetTypeE) {
551  case JetType::CaloJet:
552  writeCompoundJets<reco::CaloJet>(iEvent, iSetup);
553  break;
554  case JetType::PFJet:
555  writeCompoundJets<reco::PFJet>(iEvent, iSetup);
556  break;
557  case JetType::GenJet:
558  writeCompoundJets<reco::GenJet>(iEvent, iSetup);
559  break;
560  case JetType::BasicJet:
561  writeCompoundJets<reco::BasicJet>(iEvent, iSetup);
562  break;
563  default:
564  throw cms::Exception("InvalidInput") << "invalid jet type in CompoundJetProducer\n";
565  break;
566  };
567  } else if (writeJetsWithConst_) {
568  // Write jets and new constituents.
569  writeJetsWithConstituents<reco::PFJet>(iEvent, iSetup);
570  } else {
571  switch (jetTypeE) {
572  case JetType::CaloJet:
573  writeJets<reco::CaloJet>(iEvent, iSetup);
574  break;
575  case JetType::PFJet:
576  writeJets<reco::PFJet>(iEvent, iSetup);
577  break;
578  case JetType::GenJet:
579  writeJets<reco::GenJet>(iEvent, iSetup);
580  break;
581  case JetType::TrackJet:
582  writeJets<reco::TrackJet>(iEvent, iSetup);
583  break;
584  case JetType::PFClusterJet:
585  writeJets<reco::PFClusterJet>(iEvent, iSetup);
586  break;
587  case JetType::BasicJet:
588  writeJets<reco::BasicJet>(iEvent, iSetup);
589  break;
590  default:
591  throw cms::Exception("InvalidInput") << "invalid jet type in VirtualJetProducer\n";
592  break;
593  };
594  }
595 }
596 
597 namespace {
598  template <typename T>
599  struct Area {
600  static float get(T const&) { return 0; }
601  };
602 
603  template <>
604  struct Area<reco::CaloJet> {
605  static float get(reco::CaloJet const& jet) { return jet.getSpecific().mTowersArea; }
606  };
607 } // namespace
608 
609 template <typename T>
611  // std::cout << "writeJets " << typeid(T).name()
612  // << (doRhoFastjet_ ? " doRhoFastjet " : "")
613  // << (doAreaFastjet_ ? " doAreaFastjet " : "")
614  // << (doAreaDiskApprox_ ? " doAreaDiskApprox " : "")
615  // << std::endl;
616 
617  if (doRhoFastjet_) {
618  // declare jet collection without the two jets,
619  // for unbiased background estimation.
620  std::vector<fastjet::PseudoJet> fjexcluded_jets;
621  fjexcluded_jets = fjJets_;
622 
623  if (fjexcluded_jets.size() > 2)
624  fjexcluded_jets.resize(nExclude_);
625 
626  if (doFastJetNonUniform_) {
627  auto rhos = std::make_unique<std::vector<double>>();
628  auto sigmas = std::make_unique<std::vector<double>>();
629  int nEta = puCenters_.size();
630  rhos->reserve(nEta);
631  sigmas->reserve(nEta);
632  fastjet::ClusterSequenceAreaBase const* clusterSequenceWithArea =
633  dynamic_cast<fastjet::ClusterSequenceAreaBase const*>(fjClusterSeq_.get());
634 
635  if (clusterSequenceWithArea == nullptr) {
636  if (!fjJets_.empty()) {
637  throw cms::Exception("LogicError") << "fjClusterSeq is not initialized while inputs are present\n ";
638  }
639  } else {
640  for (int ie = 0; ie < nEta; ++ie) {
641  double eta = puCenters_[ie];
642  double etamin = eta - puWidth_;
643  double etamax = eta + puWidth_;
644  fastjet::RangeDefinition range_rho(etamin, etamax);
645  fastjet::BackgroundEstimator bkgestim(*clusterSequenceWithArea, range_rho);
646  bkgestim.set_excluded_jets(fjexcluded_jets);
647  rhos->push_back(bkgestim.rho());
648  sigmas->push_back(bkgestim.sigma());
649  }
650  }
651  iEvent.put(std::move(rhos), "rhos");
652  iEvent.put(std::move(sigmas), "sigmas");
653  } else {
654  auto rho = std::make_unique<double>(0.0);
655  auto sigma = std::make_unique<double>(0.0);
656  double mean_area = 0;
657 
658  fastjet::ClusterSequenceAreaBase const* clusterSequenceWithArea =
659  dynamic_cast<fastjet::ClusterSequenceAreaBase const*>(fjClusterSeq_.get());
660  if (clusterSequenceWithArea == nullptr) {
661  if (!fjJets_.empty()) {
662  throw cms::Exception("LogicError") << "fjClusterSeq is not initialized while inputs are present\n ";
663  }
664  } else {
665  clusterSequenceWithArea->get_median_rho_and_sigma(*fjSelector_, false, *rho, *sigma, mean_area);
666  if ((*rho < 0) || (edm::isNotFinite(*rho))) {
667  edm::LogError("BadRho") << "rho value is " << *rho << " area:" << mean_area
668  << " and n_empty_jets: " << clusterSequenceWithArea->n_empty_jets(*fjSelector_)
669  << " with range " << fjSelector_->description() << ". Setting rho to rezo.";
670  *rho = 0;
671  }
672  }
673  iEvent.put(std::move(rho), "rho");
674  iEvent.put(std::move(sigma), "sigma");
675  }
676  } // doRhoFastjet_
677 
678  // produce output jet collection
679 
680  using namespace reco;
681 
682  // allocate fjJets_.size() Ts in vector
683  auto jets = std::make_unique<std::vector<T>>(fjJets_.size());
684 
685  if (!fjJets_.empty()) {
686  // Distance between jet centers and overlap area -- for disk-based area calculation
687  using RIJ = std::pair<double, double>;
688  std::vector<RIJ> rijStorage(fjJets_.size() * (fjJets_.size() / 2));
689  RIJ* rij[fjJets_.size()];
690  unsigned int k = 0;
691  for (unsigned int ijet = 0; ijet < fjJets_.size(); ++ijet) {
692  rij[ijet] = &rijStorage[k];
693  k += ijet;
694  }
695 
696  float etaJ[fjJets_.size()], phiJ[fjJets_.size()];
697 
698  auto orParam_ = 1. / rParam_;
699  // fill jets
700  [[maybe_unused]] const CaloGeometry* pGeometry = nullptr;
701  [[maybe_unused]] const HcalTopology* pTopology = nullptr;
702  if constexpr (std::is_same_v<T, reco::CaloJet>) {
703  pGeometry = &getGeometry(iSetup);
704  pTopology = &getTopology(iSetup);
705  }
706  for (unsigned int ijet = 0; ijet < fjJets_.size(); ++ijet) {
707  auto& jet = (*jets)[ijet];
708 
709  // get the fastjet jet
710  const fastjet::PseudoJet& fjJet = fjJets_[ijet];
711  // get the constituents from fastjet
712  std::vector<fastjet::PseudoJet> const& fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
713  // convert them to CandidatePtr vector
714  std::vector<CandidatePtr> const& constituents = getConstituents(fjConstituents);
715 
716  // write the specifics to the jet (simultaneously sets 4-vector, vertex).
717  // These are overridden functions that will call the appropriate
718  // specific allocator.
719  if ((applyWeight_) && (makePFJet(jetTypeE)))
720  writeSpecific(dynamic_cast<reco::PFJet&>(jet),
721  Particle::LorentzVector(fjJet.px(), fjJet.py(), fjJet.pz(), fjJet.E()),
722  vertex_,
723  constituents,
724  &weights_);
725  else {
726  if constexpr (std::is_same_v<T, reco::CaloJet>) {
728  Particle::LorentzVector(fjJet.px(), fjJet.py(), fjJet.pz(), fjJet.E()),
729  vertex_,
730  constituents,
731  *pGeometry,
732  *pTopology);
733  } else {
735  jet, Particle::LorentzVector(fjJet.px(), fjJet.py(), fjJet.pz(), fjJet.E()), vertex_, constituents);
736  }
737  }
738  phiJ[ijet] = jet.phi();
739  etaJ[ijet] = jet.eta();
740  jet.setIsWeighted(applyWeight_);
741  }
742 
743  // calcuate the jet area
744  for (unsigned int ijet = 0; ijet < fjJets_.size(); ++ijet) {
745  // calcuate the jet area
746  double jetArea = 0.0;
747  // get the fastjet jet
748  const auto& fjJet = fjJets_[ijet];
749  if (doAreaFastjet_ && fjJet.has_area()) {
750  jetArea = fjJet.area();
751  } else if (doAreaDiskApprox_) {
752  // Here it is assumed that fjJets_ is in decreasing order of pT,
753  // which should happen in FastjetJetProducer::runAlgorithm()
754  jetArea = M_PI;
755  RIJ* distance = rij[ijet];
756  for (unsigned jJet = 0; jJet < ijet; ++jJet) {
757  distance[jJet].first = std::sqrt(reco::deltaR2(etaJ[ijet], phiJ[ijet], etaJ[jJet], phiJ[jJet])) * orParam_;
759  jetArea -= distance[jJet].second;
760  for (unsigned kJet = 0; kJet < jJet; ++kJet) {
762  distance[kJet].first,
763  rij[jJet][kJet].first,
764  distance[jJet].second,
765  distance[kJet].second,
766  rij[jJet][kJet].second);
767  } // end loop over harder jets
768  } // end loop over harder jets
769  jetArea *= (rParam_ * rParam_);
770  }
771  auto& jet = (*jets)[ijet];
772  jet.setJetArea(jetArea);
773 
774  if (doPUOffsetCorr_) {
775  jet.setPileup(subtractor_->getPileUpEnergy(ijet));
776  } else {
777  jet.setPileup(0.0);
778  }
779 
780  // std::cout << "area " << ijet << " " << jetArea << " " << Area<T>::get(jet) << std::endl;
781  // std::cout << "JetVI " << ijet << ' '<< jet.pt() << " " << jet.et() << ' '<< jet.energy() << ' '<< jet.mass() << std::endl;
782  }
783  }
784  // put the jets in the collection
785  iEvent.put(std::move(jets), jetCollInstanceName_);
786 }
787 
789 template <class T>
791  if (verbosity_ >= 1) {
792  std::cout << "<VirtualJetProducer::writeCompoundJets (moduleLabel = " << moduleLabel_ << ")>:" << std::endl;
793  }
794 
795  // get a list of output jets
796  auto jetCollection = std::make_unique<reco::BasicJetCollection>();
797  // get a list of output subjets
798  auto subjetCollection = std::make_unique<std::vector<T>>();
799 
800  // This will store the handle for the subjets after we write them
801  edm::OrphanHandle<std::vector<T>> subjetHandleAfterPut;
802  // this is the mapping of subjet to hard jet
803  std::vector<std::vector<int>> indices;
804  // this is the list of hardjet 4-momenta
805  std::vector<math::XYZTLorentzVector> p4_hardJets;
806  // this is the hardjet areas
807  std::vector<double> area_hardJets;
808 
809  // Loop over the hard jets
810  std::vector<fastjet::PseudoJet>::const_iterator it = fjJets_.begin(), iEnd = fjJets_.end(), iBegin = fjJets_.begin();
811  indices.resize(fjJets_.size());
812  for (; it != iEnd; ++it) {
813  fastjet::PseudoJet const& localJet = *it;
814  unsigned int jetIndex = it - iBegin;
815  // Get the 4-vector for the hard jet
816  p4_hardJets.push_back(math::XYZTLorentzVector(localJet.px(), localJet.py(), localJet.pz(), localJet.e()));
817  double localJetArea = 0.0;
818  if (doAreaFastjet_ && localJet.has_area()) {
819  localJetArea = localJet.area();
820  }
821  area_hardJets.push_back(localJetArea);
822 
823  // create the subjet list
824  std::vector<fastjet::PseudoJet> constituents;
825  if (it->has_pieces()) {
826  constituents = it->pieces();
827  } else if (it->has_constituents()) {
828  constituents = it->constituents();
829  }
830 
831  std::vector<fastjet::PseudoJet>::const_iterator itSubJetBegin = constituents.begin(), itSubJet = itSubJetBegin,
832  itSubJetEnd = constituents.end();
833  for (; itSubJet != itSubJetEnd; ++itSubJet) {
834  fastjet::PseudoJet const& subjet = *itSubJet;
835  if (verbosity_ >= 1) {
836  std::cout << "subjet #" << (itSubJet - itSubJetBegin) << ": Pt = " << subjet.pt() << ", eta = " << subjet.eta()
837  << ", phi = " << subjet.phi() << ", mass = " << subjet.m()
838  << " (#constituents = " << subjet.constituents().size() << ")" << std::endl;
839  std::vector<fastjet::PseudoJet> subjet_constituents = subjet.constituents();
840  int idx_constituent = 0;
841  for (std::vector<fastjet::PseudoJet>::const_iterator constituent = subjet_constituents.begin();
842  constituent != subjet_constituents.end();
843  ++constituent) {
844  if (constituent->pt() < 1.e-3)
845  continue; // CV: skip ghosts
846  std::cout << " constituent #" << idx_constituent << ": Pt = " << constituent->pt()
847  << ", eta = " << constituent->eta() << ", phi = " << constituent->phi() << ","
848  << " mass = " << constituent->m() << std::endl;
849  ++idx_constituent;
850  }
851  }
852 
853  if (verbosity_ >= 1) {
854  std::cout << "subjet #" << (itSubJet - itSubJetBegin) << ": Pt = " << subjet.pt() << ", eta = " << subjet.eta()
855  << ", phi = " << subjet.phi() << ", mass = " << subjet.m()
856  << " (#constituents = " << subjet.constituents().size() << ")" << std::endl;
857  std::vector<fastjet::PseudoJet> subjet_constituents = subjet.constituents();
858  int idx_constituent = 0;
859  for (std::vector<fastjet::PseudoJet>::const_iterator constituent = subjet_constituents.begin();
860  constituent != subjet_constituents.end();
861  ++constituent) {
862  if (constituent->pt() < 1.e-3)
863  continue; // CV: skip ghosts
864  std::cout << " constituent #" << idx_constituent << ": Pt = " << constituent->pt()
865  << ", eta = " << constituent->eta() << ", phi = " << constituent->phi() << ","
866  << " mass = " << constituent->m() << std::endl;
867  ++idx_constituent;
868  }
869  }
870 
871  math::XYZTLorentzVector p4Subjet(subjet.px(), subjet.py(), subjet.pz(), subjet.e());
872  reco::Particle::Point point(0, 0, 0);
873 
874  // This will hold ptr's to the subjets
875  std::vector<reco::CandidatePtr> subjetConstituents;
876 
877  // Get the transient subjet constituents from fastjet
878  std::vector<fastjet::PseudoJet> subjetFastjetConstituents = subjet.constituents();
879  std::vector<reco::CandidatePtr> constituents = getConstituents(subjetFastjetConstituents);
880 
881  indices[jetIndex].push_back(subjetCollection->size());
882 
883  // Add the concrete subjet type to the subjet list to write to event record
884  subjetCollection->push_back(*std::make_unique<T>());
885  auto& jet = subjetCollection->back();
886  if ((applyWeight_) && (makePFJet(jetTypeE)))
887  reco::writeSpecific(dynamic_cast<reco::PFJet&>(jet), p4Subjet, point, constituents, &weights_);
888  else if constexpr (std::is_same_v<T, reco::CaloJet>) {
890  jet, p4Subjet, point, constituents, iSetup.getData(geometry_token_), iSetup.getData(topology_token_));
891  } else {
892  reco::writeSpecific(jet, p4Subjet, point, constituents);
893  }
894  jet.setIsWeighted(applyWeight_);
895  double subjetArea = 0.0;
896  if (doAreaFastjet_ && itSubJet->has_area()) {
897  subjetArea = itSubJet->area();
898  }
899  jet.setJetArea(subjetArea);
900  }
901  }
902 
903  // put subjets into event record
904  subjetHandleAfterPut = iEvent.put(std::move(subjetCollection), jetCollInstanceName_);
905 
906  // Now create the hard jets with ptr's to the subjets as constituents
907  std::vector<math::XYZTLorentzVector>::const_iterator ip4 = p4_hardJets.begin(), ip4Begin = p4_hardJets.begin(),
908  ip4End = p4_hardJets.end();
909 
910  for (; ip4 != ip4End; ++ip4) {
911  int p4_index = ip4 - ip4Begin;
912  std::vector<int>& ind = indices[p4_index];
913  std::vector<reco::CandidatePtr> i_hardJetConstituents;
914  // Add the subjets to the hard jet
915  for (std::vector<int>::const_iterator isub = ind.begin(); isub != ind.end(); ++isub) {
916  reco::CandidatePtr candPtr(subjetHandleAfterPut, *isub, false);
917  i_hardJetConstituents.push_back(candPtr);
918  }
919  reco::Particle::Point point(0, 0, 0);
920  reco::BasicJet toput(*ip4, point, i_hardJetConstituents);
921  toput.setJetArea(area_hardJets[ip4 - ip4Begin]);
922  jetCollection->push_back(toput);
923  }
924 
925  // put hard jets into event record
926  // Store the Orphan handle for adding HTT information
928 
929  if (fromHTTTopJetProducer_) {
930  addHTTTopJetTagInfoCollection(iEvent, iSetup, oh);
931  }
932 }
933 
935 template <class T>
937  if (verbosity_ >= 1) {
938  std::cout << "<VirtualJetProducer::writeJetsWithConstituents (moduleLabel = " << moduleLabel_ << ")>:" << std::endl;
939  }
940 
941  // get a list of output jets MV: make this compatible with template
942  auto jetCollection = std::make_unique<reco::PFJetCollection>();
943 
944  // this is the mapping of jet to constituents
945  std::vector<std::vector<int>> indices;
946  // this is the list of jet 4-momenta
947  std::vector<math::XYZTLorentzVector> p4_Jets;
948  // this is the jet areas
949  std::vector<double> area_Jets;
950 
951  // get a list of output constituents
952  auto constituentCollection = std::make_unique<reco::PFCandidateCollection>();
953 
954  // This will store the handle for the constituents after we write them
955  edm::OrphanHandle<reco::PFCandidateCollection> constituentHandleAfterPut;
956 
957  // Loop over the jets and extract constituents
958  std::vector<fastjet::PseudoJet> constituentsSub;
959  std::vector<fastjet::PseudoJet>::const_iterator it = fjJets_.begin(), iEnd = fjJets_.end(), iBegin = fjJets_.begin();
960  indices.resize(fjJets_.size());
961 
962  for (; it != iEnd; ++it) {
963  fastjet::PseudoJet const& localJet = *it;
964  unsigned int jetIndex = it - iBegin;
965  // Get the 4-vector for the hard jet
966  p4_Jets.push_back(math::XYZTLorentzVector(localJet.px(), localJet.py(), localJet.pz(), localJet.e()));
967  double localJetArea = 0.0;
968  if (doAreaFastjet_ && localJet.has_area()) {
969  localJetArea = localJet.area();
970  }
971  area_Jets.push_back(localJetArea);
972 
973  // create the constituent list
974  std::vector<fastjet::PseudoJet> constituents, ghosts;
975  if (it->has_pieces())
976  constituents = it->pieces();
977  else if (it->has_constituents())
978  fastjet::SelectorIsPureGhost().sift(it->constituents(), ghosts, constituents); //filter out ghosts
979 
980  //loop over constituents of jet (can be subjets or normal constituents)
981  indices[jetIndex].reserve(constituents.size());
982  constituentsSub.reserve(constituentsSub.size() + constituents.size());
983  for (fastjet::PseudoJet const& constit : constituents) {
984  indices[jetIndex].push_back(constituentsSub.size());
985  constituentsSub.push_back(constit);
986  }
987  }
988 
989  //Loop over constituents and store in the event
990  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
991  for (fastjet::PseudoJet const& constit : constituentsSub) {
992  auto orig = inputs_[constit.user_index()];
993  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(orig->pdgId());
994  reco::PFCandidate pCand(reco::PFCandidate(orig->charge(), orig->p4(), id));
996  pVec.SetPxPyPzE(constit.px(), constit.py(), constit.pz(), constit.e());
997  pCand.setP4(pVec);
998  pCand.setSourceCandidatePtr(orig->sourceCandidatePtr(0));
999  constituentCollection->push_back(pCand);
1000  }
1001  // put constituents into event record
1002  constituentHandleAfterPut = iEvent.put(std::move(constituentCollection), jetCollInstanceName_);
1003 
1004  // Now create the jets with ptr's to the constituents
1005  std::vector<math::XYZTLorentzVector>::const_iterator ip4 = p4_Jets.begin(), ip4Begin = p4_Jets.begin(),
1006  ip4End = p4_Jets.end();
1007 
1008  for (; ip4 != ip4End; ++ip4) {
1009  int p4_index = ip4 - ip4Begin;
1010  std::vector<int>& ind = indices[p4_index];
1011  std::vector<reco::CandidatePtr> i_jetConstituents;
1012  // Add the constituents to the jet
1013  for (std::vector<int>::const_iterator iconst = ind.begin(); iconst != ind.end(); ++iconst) {
1014  reco::CandidatePtr candPtr(constituentHandleAfterPut, *iconst, false);
1015  i_jetConstituents.push_back(candPtr);
1016  }
1017  if (!i_jetConstituents.empty()) { //only keep jets which have constituents after subtraction
1018  reco::Particle::Point point(0, 0, 0);
1019  reco::PFJet jet;
1020  reco::writeSpecific(jet, *ip4, point, i_jetConstituents);
1021  jet.setJetArea(area_Jets[ip4 - ip4Begin]);
1022  jetCollection->emplace_back(jet);
1023  }
1024  }
1025 
1026  // put jets into event record
1028 }
1029 
1031  return iSetup.getData(geometry_token_);
1032 }
1033 
1035  return iSetup.getData(topology_token_);
1036 }
1037 
1038 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
1041  fillDescriptionsFromVirtualJetProducer(desc);
1042  desc.add<string>("jetCollInstanceName", "");
1043 
1044  // addDefault must be used here instead of add unless
1045  // all the classes that inherit from this class redefine
1046  // the fillDescriptions function. Otherwise, the autogenerated
1047  // cfi filenames are the same and conflict.
1048  descriptions.addDefault(desc);
1049 }
1050 
1052  desc.add<edm::InputTag>("src", edm::InputTag("particleFlow"));
1053  desc.add<edm::InputTag>("srcPVs", edm::InputTag(""));
1054  desc.add<string>("jetType", "PFJet");
1055  desc.add<string>("jetAlgorithm", "AntiKt");
1056  desc.add<double>("rParam", 0.4);
1057  desc.add<double>("inputEtMin", 0.0);
1058  desc.add<double>("inputEMin", 0.0);
1059  desc.add<double>("jetPtMin", 5.);
1060  desc.add<bool>("doPVCorrection", false);
1061  desc.add<bool>("doAreaFastjet", false);
1062  desc.add<bool>("doRhoFastjet", false);
1063  desc.add<bool>("doPUOffsetCorr", false);
1064  desc.add<double>("puPtMin", 10.);
1065  desc.add<double>("nSigmaPU", 1.0);
1066  desc.add<double>("radiusPU", 0.5);
1067  desc.add<string>("subtractorName", "");
1068  desc.add<bool>("useExplicitGhosts", false);
1069  desc.add<bool>("doAreaDiskApprox", false);
1070  desc.add<double>("voronoiRfact", -0.9);
1071  desc.add<double>("Rho_EtaMax", 4.4);
1072  desc.add<double>("Ghost_EtaMax", 5.);
1073  desc.add<int>("Active_Area_Repeats", 1);
1074  desc.add<double>("GhostArea", 0.01);
1075  desc.add<bool>("restrictInputs", false);
1076  desc.add<unsigned int>("maxInputs", 1);
1077  desc.add<bool>("writeCompound", false);
1078  desc.add<bool>("writeJetsWithConst", false);
1079  desc.add<bool>("doFastJetNonUniform", false);
1080  desc.add<bool>("useDeterministicSeed", false);
1081  desc.add<unsigned int>("minSeed", 14327);
1082  desc.add<int>("verbosity", 0);
1083  desc.add<double>("puWidth", 0.);
1084  desc.add<unsigned int>("nExclude", 0);
1085  desc.add<unsigned int>("maxBadEcalCells", 9999999);
1086  desc.add<unsigned int>("maxBadHcalCells", 9999999);
1087  desc.add<unsigned int>("maxProblematicEcalCells", 9999999);
1088  desc.add<unsigned int>("maxProblematicHcalCells", 9999999);
1089  desc.add<unsigned int>("maxRecoveredEcalCells", 9999999);
1090  desc.add<unsigned int>("maxRecoveredHcalCells", 9999999);
1091  vector<double> puCentersDefault;
1092  desc.add<vector<double>>("puCenters", puCentersDefault);
1093  desc.add<bool>("applyWeight", false);
1094  desc.add<edm::InputTag>("srcWeights", edm::InputTag(""));
1095  desc.add<double>("minimumTowersFraction", 0.);
1096 }
void writeJets(edm::Event &iEvent, edm::EventSetup const &iSetup)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Definition: helper.py:1
Jets made from CaloTowers.
Definition: CaloJet.h:27
Ptr< value_type > ptrAt(size_type i) const
HcalTopology const & getTopology(edm::EventSetup const &) const
virtual std::vector< reco::CandidatePtr > getConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents)
bool empty() const
T w() const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
virtual void inputTowers()
Base class for all types of Jets.
Definition: Jet.h:20
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, CaloGeometry const &geometry, HcalTopology const &topology)
Definition: JetSpecific.cc:32
void setSourceCandidatePtr(const PFCandidatePtr &ptr)
Definition: PFCandidate.h:120
Log< level::Error, false > LogError
virtual void copyConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents, reco::Jet *jet)
Jets made from CaloTowers.
Definition: BasicJet.h:19
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void writeJetsWithConstituents(edm::Event &iEvent, edm::EventSetup const &iSetup)
function template to write out the outputs
Jets made from PFObjects.
Definition: PFJet.h:20
const std::string names[nVars_]
static std::string const input
Definition: EdmProvDump.cc:50
size_type size() const
void writeCompoundJets(edm::Event &iEvent, edm::EventSetup const &iSetup)
function template to write out the outputs
U second(std::pair< T, U > const &p)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
virtual bool isAnomalousTower(reco::CandidatePtr input)
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:101
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
bool operator()(const fastjet::PseudoJet &t1, const fastjet::PseudoJet &t2) const
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
virtual void makeProduces(std::string s, std::string tag="")
T sqrt(T t)
Definition: SSEVec.h:19
math::XYZPoint Point
point in the space
Definition: Particle.h:25
void set_excluded_jets(const std::vector< PseudoJet > &excluded_jets)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
CaloGeometry const & getGeometry(edm::EventSetup const &) const
#define M_PI
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
double rho()
synonym for median rho [[do we have this? Both?]]
static const char *const names[]
std::shared_ptr< fastjet::JetDefinition::Plugin > PluginPtr
double sigma()
get the sigma
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
fixed size matrix
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
ParticleType translatePdgIdToType(int pdgid) const
Definition: PFCandidate.cc:231
HLT enums.
VirtualJetProducer(const edm::ParameterSet &iConfig)
static void fillDescriptionsFromVirtualJetProducer(edm::ParameterSetDescription &desc)
virtual void output(edm::Event &iEvent, edm::EventSetup const &iSetup)
#define get
Log< level::Warning, false > LogWarning
long double T
static Type byName(const std::string &name)
void setP4(const LorentzVector &p4) final
set 4-momentum
def move(src, dest)
Definition: eostools.py:511
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
math::PtEtaPhiELorentzVectorF LorentzVector
math::XYZPoint Point
point in the space
Definition: LeafCandidate.h:27
#define LogDebug(id)