CMS 3D CMS Logo

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