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 
43 #include "fastjet/SISConePlugin.hh"
44 #include "fastjet/CMSIterativeConePlugin.hh"
45 #include "fastjet/ATLASConePlugin.hh"
46 #include "fastjet/CDFMidPointPlugin.hh"
47 
48 #include <iostream>
49 #include <memory>
50 #include <algorithm>
51 #include <limits>
52 #include <cmath>
53 #include <vdt/vdtMath.h>
54 
55 using namespace std;
56 using namespace edm;
57 
58 
59 namespace reco {
60  namespace helper {
62  bool operator()( const fastjet::PseudoJet & t1, const fastjet::PseudoJet & t2 ) const {
63  return t1.perp2() > t2.perp2();
64  }
65  };
66 
67  }
68 }
69 
70 //______________________________________________________________________________
71 const char *const VirtualJetProducer::JetType::names[] = {
72  "BasicJet","GenJet","CaloJet","PFJet","TrackJet","PFClusterJet"
73 };
74 
75 
76 //______________________________________________________________________________
79 {
80  const char *const *pos = std::find(names, names + LastJetType, name);
81  if (pos == names + LastJetType) {
82  std::string errorMessage="Requested jetType not supported: "+name+"\n";
83  throw cms::Exception("Configuration",errorMessage);
84  }
85  return (Type)(pos-names);
86 }
87 
88 
90 {
91 
92 
93  if ( writeCompound_ ) {
94  produces<reco::BasicJetCollection>();
95  }
96 
97  if (makeCaloJet(jetTypeE)) {
98  produces<reco::CaloJetCollection>(tag).setBranchAlias(alias);
99  }
100  else if (makePFJet(jetTypeE)) {
101  produces<reco::PFJetCollection>(tag).setBranchAlias(alias);
102  }
103  else if (makeGenJet(jetTypeE)) {
104  produces<reco::GenJetCollection>(tag).setBranchAlias(alias);
105  }
106  else if (makeTrackJet(jetTypeE)) {
107  produces<reco::TrackJetCollection>(tag).setBranchAlias(alias);
108  }
109  else if (makePFClusterJet(jetTypeE)) {
110  produces<reco::PFClusterJetCollection>(tag).setBranchAlias(alias);
111  }
112  else if (makeBasicJet(jetTypeE)) {
113  produces<reco::BasicJetCollection>(tag).setBranchAlias(alias);
114  }
115 }
116 
118 // construction / destruction
120 
121 //______________________________________________________________________________
123 
124  moduleLabel_ = iConfig.getParameter<string> ("@module_label");
125  src_ = iConfig.getParameter<edm::InputTag>("src");
126  srcPVs_ = iConfig.getParameter<edm::InputTag>("srcPVs");
127  jetType_ = iConfig.getParameter<string> ("jetType");
128  jetAlgorithm_ = iConfig.getParameter<string> ("jetAlgorithm");
129  rParam_ = iConfig.getParameter<double> ("rParam");
130  inputEtMin_ = iConfig.getParameter<double> ("inputEtMin");
131  inputEMin_ = iConfig.getParameter<double> ("inputEMin");
132  jetPtMin_ = iConfig.getParameter<double> ("jetPtMin");
133  doPVCorrection_ = iConfig.getParameter<bool> ("doPVCorrection");
134  doAreaFastjet_ = iConfig.getParameter<bool> ("doAreaFastjet");
135  doRhoFastjet_ = iConfig.getParameter<bool> ("doRhoFastjet");
136  jetCollInstanceName_ = iConfig.getParameter<string> ("jetCollInstanceName");
137  doPUOffsetCorr_ = iConfig.getParameter<bool> ("doPUOffsetCorr");
138  puSubtractorName_ = iConfig.getParameter<string> ("subtractorName");
139  useExplicitGhosts_ = iConfig.getParameter<bool> ("useExplicitGhosts"); // use explicit ghosts in the fastjet clustering sequence?
140  doAreaDiskApprox_ = iConfig.getParameter<bool> ("doAreaDiskApprox");
141  voronoiRfact_ = iConfig.getParameter<double> ("voronoiRfact"); // Voronoi-based area calculation allows for an empirical scale factor
142  rhoEtaMax_ = iConfig.getParameter<double> ("Rho_EtaMax"); // do fasjet area / rho calcluation? => accept corresponding parameters
143  ghostEtaMax_ = iConfig.getParameter<double> ("Ghost_EtaMax");
144  activeAreaRepeats_ = iConfig.getParameter<int> ("Active_Area_Repeats");
145  ghostArea_ = iConfig.getParameter<double> ("GhostArea");
146  restrictInputs_ = iConfig.getParameter<bool> ("restrictInputs"); // restrict inputs to first "maxInputs" towers?
147  maxInputs_ = iConfig.getParameter<unsigned int>("maxInputs");
148  writeCompound_ = iConfig.getParameter<bool> ("writeCompound"); // Check to see if we are writing compound jets for substructure and jet grooming
149  doFastJetNonUniform_ = iConfig.getParameter<bool> ("doFastJetNonUniform");
150  puCenters_ = iConfig.getParameter<vector<double> >("puCenters");
151  puWidth_ = iConfig.getParameter<double> ("puWidth");
152  nExclude_ = iConfig.getParameter<unsigned int>("nExclude");
153  useDeterministicSeed_ = iConfig.getParameter<bool> ("useDeterministicSeed");
154  minSeed_ = iConfig.getParameter<unsigned int>("minSeed");
155  verbosity_ = iConfig.getParameter<int> ("verbosity");
156 
157  anomalousTowerDef_ = auto_ptr<AnomalousTower>(new AnomalousTower(iConfig));
158 
159  input_vertex_token_ = consumes<reco::VertexCollection>(srcPVs_);
160  input_candidateview_token_ = consumes<reco::CandidateView>(src_);
161  input_candidatefwdptr_token_ = consumes<vector<edm::FwdPtr<reco::PFCandidate> > >(iConfig.getParameter<edm::InputTag>("src"));
162  input_packedcandidatefwdptr_token_ = consumes<vector<edm::FwdPtr<pat::PackedCandidate> > >(iConfig.getParameter<edm::InputTag>("src"));
163  //
164  // additional parameters to think about:
165  // - overlap threshold (set to 0.75 for the time being)
166  // - p parameter for generalized kT (set to -2 for the time being)
167  // - fastjet PU subtraction parameters (not yet considered)
168  //
169  if (jetAlgorithm_=="Kt")
170  fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::kt_algorithm,rParam_));
171 
172  else if (jetAlgorithm_=="CambridgeAachen")
173  fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::cambridge_algorithm,rParam_) );
174 
175  else if (jetAlgorithm_=="AntiKt")
176  fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::antikt_algorithm,rParam_) );
177 
178  else if (jetAlgorithm_=="GeneralizedKt")
179  fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::genkt_algorithm,rParam_,-2) );
180 
181  else if (jetAlgorithm_=="SISCone") {
182 
183  fjPlugin_ = PluginPtr( new fastjet::SISConePlugin(rParam_,0.75,0,0.0,false,fastjet::SISConePlugin::SM_pttilde) );
184  fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(&*fjPlugin_) );
185 
186  } else if (jetAlgorithm_=="IterativeCone") {
187 
188  fjPlugin_ = PluginPtr(new fastjet::CMSIterativeConePlugin(rParam_,1.0));
189  fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
190 
191  } else if (jetAlgorithm_=="CDFMidPoint") {
192 
193  fjPlugin_ = PluginPtr(new fastjet::CDFMidPointPlugin(rParam_,0.75));
194  fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
195 
196  } else if (jetAlgorithm_=="ATLASCone") {
197 
198  fjPlugin_ = PluginPtr(new fastjet::ATLASConePlugin(rParam_));
199  fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
200 
201  } else {
202  throw cms::Exception("Invalid jetAlgorithm") <<"Jet algorithm for VirtualJetProducer is invalid, Abort!\n";
203  }
204 
205  jetTypeE=JetType::byName(jetType_);
206 
207  if ( doPUOffsetCorr_ ) {
208  if(puSubtractorName_.empty()){
209  LogWarning("VirtualJetProducer") << "Pile Up correction on; however, pile up type is not specified. Using default... \n";
210  subtractor_ = boost::shared_ptr<PileUpSubtractor>(new PileUpSubtractor(iConfig, consumesCollector()));
211  } else subtractor_ = boost::shared_ptr<PileUpSubtractor>(
212  PileUpSubtractorFactory::get()->create( puSubtractorName_, iConfig, consumesCollector()));
213  }
214 
215  // do approximate disk-based area calculation => warn if conflicting request
216  if (doAreaDiskApprox_ && doAreaFastjet_)
217  throw cms::Exception("Conflicting area calculations") << "Both the calculation of jet area via fastjet and via an analytical disk approximation have been requested. Please decide on one.\n";
218 
219  if ( doAreaFastjet_ || doRhoFastjet_ ) {
220 
221  if (voronoiRfact_ <= 0) {
222  fjActiveArea_ = ActiveAreaSpecPtr(new fastjet::GhostedAreaSpec(ghostEtaMax_,activeAreaRepeats_,ghostArea_));
223  fjActiveArea_->set_fj2_placement(true);
224 
225  if ( !useExplicitGhosts_ ) {
226  fjAreaDefinition_ = AreaDefinitionPtr( new fastjet::AreaDefinition(fastjet::active_area, *fjActiveArea_ ) );
227  } else {
228  fjAreaDefinition_ = AreaDefinitionPtr( new fastjet::AreaDefinition(fastjet::active_area_explicit_ghosts, *fjActiveArea_ ) );
229  }
230  }
231  fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(rhoEtaMax_) );
232  }
233 
234  if( ( doFastJetNonUniform_ ) && ( puCenters_.size() == 0 ) )
235  throw cms::Exception("doFastJetNonUniform") << "Parameter puCenters for doFastJetNonUniform is not defined." << std::endl;
236 
237  // make the "produces" statements
238  makeProduces( moduleLabel_, jetCollInstanceName_ );
239  produces<vector<double> >("rhos");
240  produces<vector<double> >("sigmas");
241  produces<double>("rho");
242  produces<double>("sigma");
243 
244 
245 }
246 
247 //______________________________________________________________________________
249 {
250 }
251 
252 
254 // implementation of member functions
256 
257 //______________________________________________________________________________
259 {
260 
261  // If requested, set the fastjet random seed to a deterministic function
262  // of the run/lumi/event.
263  // NOTE!!! The fastjet random number sequence is a global singleton.
264  // Thus, we have to create an object and get access to the global singleton
265  // in order to change it.
266  if ( useDeterministicSeed_ ) {
267  fastjet::GhostedAreaSpec gas;
268  std::vector<int> seeds(2);
269  unsigned int runNum_uint = static_cast <unsigned int> (iEvent.id().run());
270  unsigned int evNum_uint = static_cast <unsigned int> (iEvent.id().event());
271  seeds[0] = std::max(runNum_uint,minSeed_ + 3) + 3 * evNum_uint;
272  seeds[1] = std::max(runNum_uint,minSeed_ + 5) + 5 * evNum_uint;
273  gas.set_random_status(seeds);
274  }
275 
276  LogDebug("VirtualJetProducer") << "Entered produce\n";
277  //determine signal vertex2
278  vertex_=reco::Jet::Point(0,0,0);
279  if ( (makeCaloJet(jetTypeE) || makePFJet(jetTypeE)) &&doPVCorrection_) {
280  LogDebug("VirtualJetProducer") << "Adding PV info\n";
282  iEvent.getByToken(input_vertex_token_ , pvCollection);
283  if (pvCollection->size()>0) vertex_=pvCollection->begin()->position();
284  }
285 
286  // For Pileup subtraction using offset correction:
287  // set up geometry map
288  if ( doPUOffsetCorr_ ) {
289  subtractor_->setupGeometryMap(iEvent, iSetup);
290  }
291 
292  // clear data
293  LogDebug("VirtualJetProducer") << "Clear data\n";
294  fjInputs_.clear();
295  fjJets_.clear();
296  inputs_.clear();
297 
298  // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
300 
302  edm::Handle< std::vector<edm::FwdPtr<pat::PackedCandidate> > > packedinputsHandleAsFwdPtr;
303 
304  bool isView = iEvent.getByToken(input_candidateview_token_, inputsHandle);
305  if ( isView ) {
306  if ( inputsHandle->size() == 0) {
307  output( iEvent, iSetup );
308  return;
309  }
310  for (size_t i = 0; i < inputsHandle->size(); ++i) {
311  inputs_.push_back(inputsHandle->ptrAt(i));
312  }
313  } else {
314  bool isPF = iEvent.getByToken(input_candidatefwdptr_token_, pfinputsHandleAsFwdPtr);
315  if ( isPF ) {
316  if ( pfinputsHandleAsFwdPtr->size() == 0) {
317  output( iEvent, iSetup );
318  return;
319  }
320  for (size_t i = 0; i < pfinputsHandleAsFwdPtr->size(); ++i) {
321  if ( (*pfinputsHandleAsFwdPtr)[i].ptr().isAvailable() ) {
322  inputs_.push_back( (*pfinputsHandleAsFwdPtr)[i].ptr() );
323  }
324  else if ( (*pfinputsHandleAsFwdPtr)[i].backPtr().isAvailable() ) {
325  inputs_.push_back( (*pfinputsHandleAsFwdPtr)[i].backPtr() );
326  }
327  }
328  } else {
329  iEvent.getByToken(input_packedcandidatefwdptr_token_, packedinputsHandleAsFwdPtr);
330  if ( packedinputsHandleAsFwdPtr->size() == 0) {
331  output( iEvent, iSetup );
332  return;
333  }
334  for (size_t i = 0; i < packedinputsHandleAsFwdPtr->size(); ++i) {
335  if ( (*packedinputsHandleAsFwdPtr)[i].ptr().isAvailable() ) {
336  inputs_.push_back( (*packedinputsHandleAsFwdPtr)[i].ptr() );
337  }
338  else if ( (*packedinputsHandleAsFwdPtr)[i].backPtr().isAvailable() ) {
339  inputs_.push_back( (*packedinputsHandleAsFwdPtr)[i].backPtr() );
340  }
341  }
342  }
343  }
344  LogDebug("VirtualJetProducer") << "Got inputs\n";
345 
346  // Convert candidates to fastjet::PseudoJets.
347  // Also correct to Primary Vertex. Will modify fjInputs_
348  // and use inputs_
349  fjInputs_.reserve(inputs_.size());
350  inputTowers();
351  LogDebug("VirtualJetProducer") << "Inputted towers\n";
352 
353  // For Pileup subtraction using offset correction:
354  // Subtract pedestal.
355  if ( doPUOffsetCorr_ ) {
356  subtractor_->setDefinition(fjJetDefinition_);
357  subtractor_->reset(inputs_,fjInputs_,fjJets_);
358  subtractor_->calculatePedestal(fjInputs_);
359  subtractor_->subtractPedestal(fjInputs_);
360  LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
361  }
362  // Run algorithm. Will modify fjJets_ and allocate fjClusterSeq_.
363  // This will use fjInputs_
364  runAlgorithm( iEvent, iSetup );
365 
366  // if ( doPUOffsetCorr_ ) {
367  // subtractor_->setAlgorithm(fjClusterSeq_);
368  // }
369 
370  LogDebug("VirtualJetProducer") << "Ran algorithm\n";
371  // For Pileup subtraction using offset correction:
372  // Now we find jets and need to recalculate their energy,
373  // mark towers participated in jet,
374  // remove occupied towers from the list and recalculate mean and sigma
375  // put the initial towers collection to the jet,
376  // and subtract from initial towers in jet recalculated mean and sigma of towers
377  if ( doPUOffsetCorr_ ) {
378  LogDebug("VirtualJetProducer") << "Do PUOffsetCorr\n";
379  vector<fastjet::PseudoJet> orphanInput;
380  subtractor_->calculateOrphanInput(orphanInput);
381  subtractor_->calculatePedestal(orphanInput);
382  subtractor_->offsetCorrectJets();
383  }
384  // Write the output jets.
385  // This will (by default) call the member function template
386  // "writeJets", but can be overridden.
387  // this will use inputs_
388  output( iEvent, iSetup );
389  LogDebug("VirtualJetProducer") << "Wrote jets\n";
390 
391  // Clear the work vectors so that memory is free for other modules.
392  // Use the trick of swapping with an empty vector so that the memory
393  // is actually given back rather than silently kept.
394  decltype(fjInputs_)().swap(fjInputs_);
395  decltype(fjJets_)().swap(fjJets_);
396  decltype(inputs_)().swap(inputs_);
397 
398  return;
399 }
400 
401 //______________________________________________________________________________
402 
404 {
405  auto inBegin = inputs_.begin(),
406  inEnd = inputs_.end(), i = inBegin;
407  for (; i != inEnd; ++i ) {
408  auto const & input = **i;
409  // std::cout << "CaloTowerVI jets " << input->pt() << " " << input->et() << ' '<< input->energy() << ' ' << (isAnomalousTower(input) ? " bad" : " ok") << std::endl;
410  if (edm::isNotFinite(input.pt())) continue;
411  if (input.et() <inputEtMin_) continue;
412  if (input.energy()<inputEMin_) continue;
413  if (isAnomalousTower(*i)) continue;
414  // Change by SRR : this is no longer an error nor warning, this can happen with PU mitigation algos.
415  // Also switch to something more numerically safe. (VI: 10^-42GeV????)
416  if (input.pt() < 100 * std::numeric_limits<double>::epsilon() ) {
417  continue;
418  }
419  if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
420  const CaloTower & tower = dynamic_cast<const CaloTower &>(input);
421  auto const & ct = tower.p4(vertex_); // very expensive as computed in eta/phi
422  fjInputs_.emplace_back(ct.px(),ct.py(),ct.pz(),ct.energy());
423  //std::cout << "tower:" << *tower << '\n';
424  }
425  else {
426  /*
427  if(makePFJet(jetTypeE)) {
428  reco::PFCandidate& pfc = (reco::PFCandidate&)input;
429  std::cout << "PF cand:" << pfc << '\n';
430  }
431  */
432  fjInputs_.emplace_back(input.px(),input.py(),input.pz(),
433  input.energy());
434  }
435  fjInputs_.back().set_user_index(i - inBegin);
436  }
437 
438  if ( restrictInputs_ && fjInputs_.size() > maxInputs_ ) {
440  std::sort(fjInputs_.begin(), fjInputs_.end(), pTComparator);
441  fjInputs_.resize(maxInputs_);
442  edm::LogWarning("JetRecoTooManyEntries") << "Too many inputs in the event, limiting to first " << maxInputs_ << ". Output is suspect.";
443  }
444 }
445 
446 //______________________________________________________________________________
448 {
449  if (!makeCaloJet(jetTypeE))
450  return false;
451  else
452  return (*anomalousTowerDef_)(*input);
453 }
454 
455 //------------------------------------------------------------------------------
456 // This is pure virtual.
457 //______________________________________________________________________________
458 // void VirtualJetProducer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup,
459 // std::vector<edm::Ptr<reco::Candidate> > const & inputs_);
460 
461 //______________________________________________________________________________
462 void VirtualJetProducer::copyConstituents(const vector<fastjet::PseudoJet>& fjConstituents,
463  reco::Jet* jet)
464 {
465  for (unsigned int i=0;i<fjConstituents.size();++i) {
466  int index = fjConstituents[i].user_index();
467  if ( index >= 0 && static_cast<unsigned int>(index) < inputs_.size() )
468  jet->addDaughter(inputs_[index]);
469  }
470 }
471 
472 
473 //______________________________________________________________________________
474 vector<reco::CandidatePtr>
475 VirtualJetProducer::getConstituents(const vector<fastjet::PseudoJet>&fjConstituents)
476 {
477  vector<reco::CandidatePtr> result; result.reserve(fjConstituents.size()/2);
478  for (unsigned int i=0;i<fjConstituents.size();i++) {
479  auto index = fjConstituents[i].user_index();
480  if ( index >= 0 && static_cast<unsigned int>(index) < inputs_.size() ) {
481  result.emplace_back(inputs_[index]);
482  }
483  }
484  return result;
485 }
486 
487 
488 //_____________________________________________________________________________
489 
491 {
492  // Write jets and constitutents. Will use fjJets_, inputs_
493  // and fjClusterSeq_
494 
495  if ( !writeCompound_ ) {
496  switch( jetTypeE ) {
497  case JetType::CaloJet :
498  writeJets<reco::CaloJet>( iEvent, iSetup);
499  break;
500  case JetType::PFJet :
501  writeJets<reco::PFJet>( iEvent, iSetup);
502  break;
503  case JetType::GenJet :
504  writeJets<reco::GenJet>( iEvent, iSetup);
505  break;
506  case JetType::TrackJet :
507  writeJets<reco::TrackJet>( iEvent, iSetup);
508  break;
509  case JetType::PFClusterJet :
510  writeJets<reco::PFClusterJet>( iEvent, iSetup);
511  break;
512  case JetType::BasicJet :
513  writeJets<reco::BasicJet>( iEvent, iSetup);
514  break;
515  default:
516  throw cms::Exception("InvalidInput") << "invalid jet type in VirtualJetProducer\n";
517  break;
518  };
519  } else {
520  // Write jets and constitutents.
521  switch( jetTypeE ) {
522  case JetType::CaloJet :
523  writeCompoundJets<reco::CaloJet>( iEvent, iSetup );
524  break;
525  case JetType::PFJet :
526  writeCompoundJets<reco::PFJet>( iEvent, iSetup );
527  break;
528  case JetType::GenJet :
529  writeCompoundJets<reco::GenJet>( iEvent, iSetup );
530  break;
531  case JetType::BasicJet :
532  writeCompoundJets<reco::BasicJet>( iEvent, iSetup );
533  break;
534  default:
535  throw cms::Exception("InvalidInput") << "invalid jet type in CompoundJetProducer\n";
536  break;
537  };
538  }
539 }
540 
541 namespace {
542 template< typename T >
543 struct Area { static float get(T const &) {return 0;}};
544 
545 template<>
546 struct Area<reco::CaloJet>{ static float get(reco::CaloJet const & jet) {
547  return jet.getSpecific().mTowersArea;
548 }
549 };
550 }
551 
552 template< typename T >
554 {
555  // std::cout << "writeJets " << typeid(T).name()
556  // << (doRhoFastjet_ ? " doRhoFastjet " : "")
557  // << (doAreaFastjet_ ? " doAreaFastjet " : "")
558  // << (doAreaDiskApprox_ ? " doAreaDiskApprox " : "")
559  // << std::endl;
560 
561  if (doRhoFastjet_) {
562  // declare jet collection without the two jets,
563  // for unbiased background estimation.
564  std::vector<fastjet::PseudoJet> fjexcluded_jets;
565  fjexcluded_jets=fjJets_;
566 
567  if(fjexcluded_jets.size()>2) fjexcluded_jets.resize(nExclude_);
568 
569  if(doFastJetNonUniform_){
570  auto rhos = std::make_unique<std::vector<double>>();
571  auto sigmas = std::make_unique<std::vector<double>>();
572  int nEta = puCenters_.size();
573  rhos->reserve(nEta);
574  sigmas->reserve(nEta);
575  fastjet::ClusterSequenceAreaBase const* clusterSequenceWithArea =
576  dynamic_cast<fastjet::ClusterSequenceAreaBase const *> ( &*fjClusterSeq_ );
577 
578  if (clusterSequenceWithArea ==nullptr ){
579  if (fjJets_.size() > 0) {
580  throw cms::Exception("LogicError")<<"fjClusterSeq is not initialized while inputs are present\n ";
581  }
582  } else {
583  for(int ie = 0; ie < nEta; ++ie){
584  double eta = puCenters_[ie];
585  double etamin=eta-puWidth_;
586  double etamax=eta+puWidth_;
587  fastjet::RangeDefinition range_rho(etamin,etamax);
588  fastjet::BackgroundEstimator bkgestim(*clusterSequenceWithArea,range_rho);
589  bkgestim.set_excluded_jets(fjexcluded_jets);
590  rhos->push_back(bkgestim.rho());
591  sigmas->push_back(bkgestim.sigma());
592  }
593  }
594  iEvent.put(std::move(rhos),"rhos");
595  iEvent.put(std::move(sigmas),"sigmas");
596  }else{
597  auto rho = std::make_unique<double>(0.0);
598  auto sigma = std::make_unique<double>(0.0);
599  double mean_area = 0;
600 
601  fastjet::ClusterSequenceAreaBase const* clusterSequenceWithArea =
602  dynamic_cast<fastjet::ClusterSequenceAreaBase const *> ( &*fjClusterSeq_ );
603  /*
604  const double nemptyjets = clusterSequenceWithArea->n_empty_jets(*fjRangeDef_);
605  if(( nemptyjets < -15 ) || ( nemptyjets > fjRangeDef_->area()+ 15)) {
606  edm::LogWarning("StrangeNEmtpyJets") << "n_empty_jets is : " << clusterSequenceWithArea->n_empty_jets(*fjRangeDef_) << " with range " << fjRangeDef_->description() << ".";
607  }
608  */
609  if (clusterSequenceWithArea ==nullptr ){
610  if (fjJets_.size() > 0) {
611  throw cms::Exception("LogicError")<<"fjClusterSeq is not initialized while inputs are present\n ";
612  }
613  } else {
614  clusterSequenceWithArea->get_median_rho_and_sigma(*fjRangeDef_,false,*rho,*sigma,mean_area);
615  if((*rho < 0)|| (edm::isNotFinite(*rho))) {
616  edm::LogError("BadRho") << "rho value is " << *rho << " area:" << mean_area << " and n_empty_jets: " << clusterSequenceWithArea->n_empty_jets(*fjRangeDef_) << " with range " << fjRangeDef_->description()
617  <<". Setting rho to rezo.";
618  *rho = 0;
619  }
620  }
621  iEvent.put(std::move(rho),"rho");
622  iEvent.put(std::move(sigma),"sigma");
623  }
624  } // doRhoFastjet_
625 
626  // produce output jet collection
627 
628  using namespace reco;
629 
630  // allocate fjJets_.size() Ts in vector
631  auto jets = std::make_unique<std::vector<T>>(fjJets_.size());
632 
633  // Distance between jet centers and overlap area -- for disk-based area calculation
634  using RIJ = std::pair<double,double>;
635  std::vector<RIJ> rijStorage(fjJets_.size()*(fjJets_.size()/2));
636  RIJ * rij[fjJets_.size()];
637  unsigned int k=0;
638  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
639  rij[ijet] = &rijStorage[k]; k+=ijet;
640  }
641 
642  float etaJ[fjJets_.size()], phiJ[fjJets_.size()];
643 
644  auto orParam_ = 1./rParam_;
645  // fill jets
646  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
647  auto & jet = (*jets)[ijet];
648  // get the fastjet jet
649  const fastjet::PseudoJet& fjJet = fjJets_[ijet];
650  // get the constituents from fastjet
651  std::vector<fastjet::PseudoJet> const & fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
652  // convert them to CandidatePtr vector
653  std::vector<CandidatePtr> const & constituents = getConstituents(fjConstituents);
654 
655  // write the specifics to the jet (simultaneously sets 4-vector, vertex).
656  // These are overridden functions that will call the appropriate
657  // specific allocator.
659  Particle::LorentzVector(fjJet.px(),
660  fjJet.py(),
661  fjJet.pz(),
662  fjJet.E()),
663  vertex_,
664  constituents, iSetup);
665  phiJ[ijet] = jet.phi();
666  etaJ[ijet] = jet.eta();
667  }
668 
669  // calcuate the jet area
670  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
671  // calcuate the jet area
672  double jetArea=0.0;
673  // get the fastjet jet
674  const auto & fjJet = fjJets_[ijet];
675  if ( doAreaFastjet_ && fjJet.has_area() ) {
676  jetArea = fjJet.area();
677  }
678  else if ( doAreaDiskApprox_ ) {
679  // Here it is assumed that fjJets_ is in decreasing order of pT,
680  // which should happen in FastjetJetProducer::runAlgorithm()
681  jetArea = M_PI;
682  RIJ * distance = rij[ijet];
683  for (unsigned jJet = 0; jJet < ijet; ++jJet) {
684  distance[jJet].first = std::sqrt(reco::deltaR2(etaJ[ijet],phiJ[ijet], etaJ[jJet],phiJ[jJet]))*orParam_;
685  distance[jJet].second = reco::helper::VirtualJetProducerHelper::intersection(distance[jJet].first);
686  jetArea -=distance[jJet].second;
687  for (unsigned kJet = 0; kJet < jJet; ++kJet) {
688  jetArea += reco::helper::VirtualJetProducerHelper::intersection(distance[jJet].first, distance[kJet].first, rij[jJet][kJet].first,
689  distance[jJet].second, distance[kJet].second, rij[jJet][kJet].second);
690  } // end loop over harder jets
691  } // end loop over harder jets
692  jetArea *= (rParam_*rParam_);
693  }
694  auto & jet = (*jets)[ijet];
695  jet.setJetArea (jetArea);
696 
697  if(doPUOffsetCorr_){
698  jet.setPileup(subtractor_->getPileUpEnergy(ijet));
699  }else{
700  jet.setPileup (0.0);
701  }
702 
703  // std::cout << "area " << ijet << " " << jetArea << " " << Area<T>::get(jet) << std::endl;
704  // std::cout << "JetVI " << ijet << ' '<< jet.pt() << " " << jet.et() << ' '<< jet.energy() << ' '<< jet.mass() << std::endl;
705 
706  }
707  // put the jets in the collection
708  iEvent.put(std::move(jets),jetCollInstanceName_);
709 }
710 
712 template< class T>
714 {
715  if ( verbosity_ >= 1 ) {
716  std::cout << "<VirtualJetProducer::writeCompoundJets (moduleLabel = " << moduleLabel_ << ")>:" << std::endl;
717  }
718 
719  // get a list of output jets
720  auto jetCollection = std::make_unique<reco::BasicJetCollection>();
721  // get a list of output subjets
722  auto subjetCollection = std::make_unique<std::vector<T>>();
723 
724  // This will store the handle for the subjets after we write them
725  edm::OrphanHandle< std::vector<T> > subjetHandleAfterPut;
726  // this is the mapping of subjet to hard jet
727  std::vector< std::vector<int> > indices;
728  // this is the list of hardjet 4-momenta
729  std::vector<math::XYZTLorentzVector> p4_hardJets;
730  // this is the hardjet areas
731  std::vector<double> area_hardJets;
732 
733  // Loop over the hard jets
734  std::vector<fastjet::PseudoJet>::const_iterator it = fjJets_.begin(),
735  iEnd = fjJets_.end(),
736  iBegin = fjJets_.begin();
737  indices.resize( fjJets_.size() );
738  for ( ; it != iEnd; ++it ) {
739  fastjet::PseudoJet const & localJet = *it;
740  unsigned int jetIndex = it - iBegin;
741  // Get the 4-vector for the hard jet
742  p4_hardJets.push_back( math::XYZTLorentzVector(localJet.px(), localJet.py(), localJet.pz(), localJet.e() ));
743  double localJetArea = 0.0;
744  if ( doAreaFastjet_ && localJet.has_area() ) {
745  localJetArea = localJet.area();
746  }
747  area_hardJets.push_back( localJetArea );
748 
749  // create the subjet list
750  std::vector<fastjet::PseudoJet> constituents;
751  if ( it->has_pieces() ) {
752  constituents = it->pieces();
753  } else if ( it->has_constituents() ) {
754  constituents = it->constituents();
755  }
756 
757  std::vector<fastjet::PseudoJet>::const_iterator itSubJetBegin = constituents.begin(),
758  itSubJet = itSubJetBegin, itSubJetEnd = constituents.end();
759  for (; itSubJet != itSubJetEnd; ++itSubJet ){
760 
761  fastjet::PseudoJet const & subjet = *itSubJet;
762  if ( verbosity_ >= 1 ) {
763  std::cout << "subjet #" << (itSubJet - itSubJetBegin) << ": Pt = " << subjet.pt() << ", eta = " << subjet.eta() << ", phi = " << subjet.phi() << ", mass = " << subjet.m()
764  << " (#constituents = " << subjet.constituents().size() << ")" << std::endl;
765  std::vector<fastjet::PseudoJet> subjet_constituents = subjet.constituents();
766  int idx_constituent = 0;
767  for ( std::vector<fastjet::PseudoJet>::const_iterator constituent = subjet_constituents.begin();
768  constituent != subjet_constituents.end(); ++constituent ) {
769  if ( constituent->pt() < 1.e-3 ) continue; // CV: skip ghosts
770  std::cout << " constituent #" << idx_constituent << ": Pt = " << constituent->pt() << ", eta = " << constituent->eta() << ", phi = " << constituent->phi() << ","
771  << " mass = " << constituent->m() << std::endl;
772  ++idx_constituent;
773  }
774  }
775 
776  if ( verbosity_ >= 1 ) {
777  std::cout << "subjet #" << (itSubJet - itSubJetBegin) << ": Pt = " << subjet.pt() << ", eta = " << subjet.eta() << ", phi = " << subjet.phi() << ", mass = " << subjet.m()
778  << " (#constituents = " << subjet.constituents().size() << ")" << std::endl;
779  std::vector<fastjet::PseudoJet> subjet_constituents = subjet.constituents();
780  int idx_constituent = 0;
781  for ( std::vector<fastjet::PseudoJet>::const_iterator constituent = subjet_constituents.begin();
782  constituent != subjet_constituents.end(); ++constituent ) {
783  if ( constituent->pt() < 1.e-3 ) continue; // CV: skip ghosts
784  std::cout << " constituent #" << idx_constituent << ": Pt = " << constituent->pt() << ", eta = " << constituent->eta() << ", phi = " << constituent->phi() << ","
785  << " mass = " << constituent->m() << std::endl;
786  ++idx_constituent;
787  }
788  }
789 
790  math::XYZTLorentzVector p4Subjet(subjet.px(), subjet.py(), subjet.pz(), subjet.e() );
792 
793  // This will hold ptr's to the subjets
794  std::vector<reco::CandidatePtr> subjetConstituents;
795 
796  // Get the transient subjet constituents from fastjet
797  std::vector<fastjet::PseudoJet> subjetFastjetConstituents = subjet.constituents();
798  std::vector<reco::CandidatePtr> constituents =
799  getConstituents(subjetFastjetConstituents );
800 
801  indices[jetIndex].push_back( subjetCollection->size() );
802 
803  // Add the concrete subjet type to the subjet list to write to event record
804  T jet;
805  reco::writeSpecific( jet, p4Subjet, point, constituents, iSetup);
806  double subjetArea = 0.0;
807  if ( doAreaFastjet_ && itSubJet->has_area() ){
808  subjetArea = itSubJet->area();
809  }
810  jet.setJetArea( subjetArea );
811  subjetCollection->push_back( jet );
812  }
813  }
814 
815  // put subjets into event record
816  subjetHandleAfterPut = iEvent.put(std::move(subjetCollection), jetCollInstanceName_);
817 
818  // Now create the hard jets with ptr's to the subjets as constituents
819  std::vector<math::XYZTLorentzVector>::const_iterator ip4 = p4_hardJets.begin(),
820  ip4Begin = p4_hardJets.begin(),
821  ip4End = p4_hardJets.end();
822 
823  for ( ; ip4 != ip4End; ++ip4 ) {
824  int p4_index = ip4 - ip4Begin;
825  std::vector<int> & ind = indices[p4_index];
826  std::vector<reco::CandidatePtr> i_hardJetConstituents;
827  // Add the subjets to the hard jet
828  for( std::vector<int>::const_iterator isub = ind.begin();
829  isub != ind.end(); ++isub ) {
830  reco::CandidatePtr candPtr( subjetHandleAfterPut, *isub, false );
831  i_hardJetConstituents.push_back( candPtr );
832  }
834  reco::BasicJet toput( *ip4, point, i_hardJetConstituents);
835  toput.setJetArea( area_hardJets[ip4 - ip4Begin] );
836  jetCollection->push_back( toput );
837  }
838 
839  // put hard jets into event record
840  // Store the Orphan handle for adding HTT information
841  edm::OrphanHandle<reco::BasicJetCollection> oh = iEvent.put(std::move(jetCollection));
842 
843  if (fromHTTTopJetProducer_){
844  addHTTTopJetTagInfoCollection( iEvent, iSetup, oh);
845  }
846 
847 }
848 
849 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
851 
853  fillDescriptionsFromVirtualJetProducer(desc);
854  desc.add<string>("jetCollInstanceName", "" );
855  descriptions.add("VirtualJetProducer",desc);
856 }
857 
859 {
860  desc.add<edm::InputTag>("src", edm::InputTag("particleFlow") );
861  desc.add<edm::InputTag>("srcPVs", edm::InputTag("") );
862  desc.add<string>("jetType", "PFJet" );
863  desc.add<string>("jetAlgorithm", "AntiKt" );
864  desc.add<double>("rParam", 0.4 );
865  desc.add<double>("inputEtMin", 0.0 );
866  desc.add<double>("inputEMin", 0.0 );
867  desc.add<double>("jetPtMin", 5. );
868  desc.add<bool> ("doPVCorrection", false );
869  desc.add<bool> ("doAreaFastjet", false );
870  desc.add<bool> ("doRhoFastjet", false );
871  desc.add<bool> ("doPUOffsetCorr", false );
872  desc.add<double>("puPtMin", 10.);
873  desc.add<double>("nSigmaPU", 1.0 );
874  desc.add<double>("radiusPU", 0.5 );
875  desc.add<string>("subtractorName", "" );
876  desc.add<bool> ("useExplicitGhosts", false );
877  desc.add<bool> ("doAreaDiskApprox", false );
878  desc.add<double>("voronoiRfact", -0.9 );
879  desc.add<double>("Rho_EtaMax", 4.4 );
880  desc.add<double>("Ghost_EtaMax", 5. );
881  desc.add<int> ("Active_Area_Repeats", 1 );
882  desc.add<double>("GhostArea", 0.01 );
883  desc.add<bool> ("restrictInputs", false );
884  desc.add<unsigned int> ("maxInputs", 1 );
885  desc.add<bool> ("writeCompound", false );
886  desc.add<bool> ("doFastJetNonUniform", false );
887  desc.add<bool> ("useDeterministicSeed",false );
888  desc.add<unsigned int> ("minSeed", 14327 );
889  desc.add<int> ("verbosity", 0 );
890  desc.add<double>("puWidth", 0. );
891  desc.add<unsigned int>("nExclude", 0 );
892  desc.add<unsigned int>("maxBadEcalCells", 9999999 );
893  desc.add<unsigned int>("maxBadHcalCells", 9999999 );
894  desc.add<unsigned int>("maxProblematicEcalCells", 9999999 );
895  desc.add<unsigned int>("maxProblematicHcalCells", 9999999 );
896  desc.add<unsigned int>("maxRecoveredEcalCells", 9999999 );
897  desc.add<unsigned int>("maxRecoveredHcalCells", 9999999 );
898  vector<double> puCentersDefault;
899  desc.add<vector<double>>("puCenters", puCentersDefault);
900 }
901 
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:39
boost::shared_ptr< fastjet::AreaDefinition > AreaDefinitionPtr
void writeJets(edm::Event &iEvent, edm::EventSetup const &iSetup)
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
boost::shared_ptr< fastjet::JetDefinition > JetDefPtr
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
Definition: helper.py:1
Jets made from CaloTowers.
Definition: CaloJet.h:29
static const HistoName names[]
math::PtEtaPhiMLorentzVector p4(double vtxZ) const
Definition: CaloTower.cc:129
virtual std::vector< reco::CandidatePtr > getConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
Ptr< value_type > ptrAt(size_type i) const
virtual void inputTowers()
Base class for all types of Jets.
Definition: Jet.h:20
boost::shared_ptr< fastjet::JetDefinition::Plugin > PluginPtr
size_type size() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual void copyConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents, reco::Jet *jet)
Jets made from CaloTowers.
Definition: BasicJet.h:20
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
static std::string const input
Definition: EdmProvDump.cc:44
void writeCompoundJets(edm::Event &iEvent, edm::EventSetup const &iSetup)
function template to write out the outputs
bool operator()(const fastjet::PseudoJet &t1, const fastjet::PseudoJet &t2) const
boost::shared_ptr< fastjet::RangeDefinition > RangeDefPtr
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:103
int iEvent
Definition: GenABIO.cc:230
bool isNotFinite(T x)
Definition: isFinite.h:10
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
virtual void makeProduces(std::string s, std::string tag="")
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
math::XYZPoint Point
point in the space
Definition: Particle.h:25
void set_excluded_jets(const std::vector< PseudoJet > &excluded_jets)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define M_PI
int k[5][pyjets_maxn]
double rho()
synonym for median rho [[do we have this? Both?]]
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static const char *const names[]
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
double sigma()
get the sigma
edm::EventID id() const
Definition: EventBase.h:60
fixed size matrix
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
HLT enums.
VirtualJetProducer(const edm::ParameterSet &iConfig)
boost::shared_ptr< fastjet::GhostedAreaSpec > ActiveAreaSpecPtr
void addDaughter(const CandidatePtr &)
add a daughter via a reference
static void fillDescriptionsFromVirtualJetProducer(edm::ParameterSetDescription &desc)
virtual void output(edm::Event &iEvent, edm::EventSetup const &iSetup)
long double T
static Type byName(const std::string &name)
def move(src, dest)
Definition: eostools.py:510
T get(const Candidate &c)
Definition: component.h:55
*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
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