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 
13 
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 *VirtualJetProducer::JetType::names[] = {
65  "BasicJet","GenJet","CaloJet","PFJet","TrackJet","PFClusterJet"
66 };
67 
68 
69 //______________________________________________________________________________
72 {
73  const char **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_ ) {
187  throw cms::Exception("InvalidInput") << "Can only offset correct jets of type CaloJet or BasicJet";
188  }
189 
190  if(iConfig.exists("subtractorName")) puSubtractorName_ = iConfig.getParameter<string> ("subtractorName");
192 
193  if(puSubtractorName_.empty()){
194  edm::LogWarning("VirtualJetProducer") << "Pile Up correction on; however, pile up type is not specified. Using default... \n";
195  subtractor_ = boost::shared_ptr<PileUpSubtractor>(new PileUpSubtractor(iConfig));
196  }else{
197  subtractor_ = boost::shared_ptr<PileUpSubtractor>(PileUpSubtractorFactory::get()->create( puSubtractorName_, iConfig));
198  }
199  }
200 
201  // use explicit ghosts in the fastjet clustering sequence?
202  if ( iConfig.exists("useExplicitGhosts") ) {
203  useExplicitGhosts_ = iConfig.getParameter<bool>("useExplicitGhosts");
204  }
205 
206  // do approximate disk-based area calculation => warn if conflicting request
207  if (iConfig.exists("doAreaDiskApprox")) {
208  doAreaDiskApprox_ = iConfig.getParameter<bool>("doAreaDiskApprox");
210  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";
211 
212  }
213  // turn off jet collection output for speed
214  // Voronoi-based area calculation allows for an empirical scale factor
215  if (iConfig.exists("voronoiRfact"))
216  voronoiRfact_ = iConfig.getParameter<double>("voronoiRfact");
217 
218 
219  // do fasjet area / rho calcluation? => accept corresponding parameters
220  if ( doAreaFastjet_ || doRhoFastjet_ ) {
221  // Eta range of jets to be considered for Rho calculation
222  // Should be at most (jet acceptance - jet radius)
223  double rhoEtaMax=iConfig.getParameter<double>("Rho_EtaMax");
224  // default Ghost_EtaMax should be 5
225  double ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
226  // default Active_Area_Repeats 1
227  int activeAreaRepeats = iConfig.getParameter<int> ("Active_Area_Repeats");
228  // default GhostArea 0.01
229  double ghostArea = iConfig.getParameter<double> ("GhostArea");
230  if (voronoiRfact_ <= 0) {
231  fjActiveArea_ = ActiveAreaSpecPtr(new fastjet::GhostedAreaSpec(ghostEtaMax,activeAreaRepeats,ghostArea));
232  fjActiveArea_->set_fj2_placement(true);
233  if ( ! useExplicitGhosts_ ) {
234  fjAreaDefinition_ = AreaDefinitionPtr( new fastjet::AreaDefinition(fastjet::active_area, *fjActiveArea_ ) );
235  } else {
236  fjAreaDefinition_ = AreaDefinitionPtr( new fastjet::AreaDefinition(fastjet::active_area_explicit_ghosts, *fjActiveArea_ ) );
237  }
238  }
239  fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(rhoEtaMax) );
240  }
241 
242  // restrict inputs to first "maxInputs" towers?
243  if ( iConfig.exists("restrictInputs") ) {
244  restrictInputs_ = iConfig.getParameter<bool>("restrictInputs");
245  maxInputs_ = iConfig.getParameter<unsigned int>("maxInputs");
246  }
247 
248 
249  string alias=iConfig.getUntrackedParameter<string>("alias",moduleLabel_);
250 
251 
252  // Check to see if we are writing compound jets for substructure
253  // and jet grooming
254  if ( iConfig.exists("writeCompound") ) {
255  writeCompound_ = iConfig.getParameter<bool>("writeCompound");
256  }
257 
258  // make the "produces" statements
260 
261  doFastJetNonUniform_ = false;
262  if(iConfig.exists("doFastJetNonUniform")) doFastJetNonUniform_ = iConfig.getParameter<bool> ("doFastJetNonUniform");
264  puCenters_ = iConfig.getParameter<std::vector<double> >("puCenters");
265  puWidth_ = iConfig.getParameter<double>("puWidth");
266  nExclude_ = iConfig.getParameter<unsigned int>("nExclude");
267  }
268 
269  useDeterministicSeed_ = false;
270  minSeed_ = 0;
271  if ( iConfig.exists("useDeterministicSeed") ) {
272  useDeterministicSeed_ = iConfig.getParameter<bool>("useDeterministicSeed");
273  minSeed_ = iConfig.getParameter<unsigned int>("minSeed");
274  }
275 
276  if ( iConfig.exists("verbosity" ) ) {
277  verbosity_ = iConfig.getParameter<int>("verbosity");
278  }
279 
280  produces<std::vector<double> >("rhos");
281  produces<std::vector<double> >("sigmas");
282  produces<double>("rho");
283  produces<double>("sigma");
284 
285 
286 }
287 
288 //______________________________________________________________________________
290 {
291 }
292 
293 
295 // implementation of member functions
297 
298 //______________________________________________________________________________
300 {
301 
302  // If requested, set the fastjet random seed to a deterministic function
303  // of the run/lumi/event.
304  // NOTE!!! The fastjet random number sequence is a global singleton.
305  // Thus, we have to create an object and get access to the global singleton
306  // in order to change it.
307  if ( useDeterministicSeed_ ) {
308  fastjet::GhostedAreaSpec gas;
309  std::vector<int> seeds(2);
310  seeds[0] = std::max(iEvent.id().run(),minSeed_ + 3) + 3 * iEvent.id().event();
311  seeds[1] = std::max(iEvent.id().run(),minSeed_ + 5) + 5 * iEvent.id().event();
312  gas.set_random_status(seeds);
313  }
314 
315  LogDebug("VirtualJetProducer") << "Entered produce\n";
316  //determine signal vertex2
317  vertex_=reco::Jet::Point(0,0,0);
319  LogDebug("VirtualJetProducer") << "Adding PV info\n";
321  iEvent.getByLabel(srcPVs_,pvCollection);
322  if (pvCollection->size()>0) vertex_=pvCollection->begin()->position();
323  }
324 
325  // For Pileup subtraction using offset correction:
326  // set up geometry map
327  if ( doPUOffsetCorr_ ) {
328  subtractor_->setupGeometryMap(iEvent, iSetup);
329  }
330 
331  // clear data
332  LogDebug("VirtualJetProducer") << "Clear data\n";
333  fjInputs_.clear();
334  fjJets_.clear();
335  inputs_.clear();
336 
337  // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
339 
341 
342  bool isView = iEvent.getByLabel(src_,inputsHandle);
343  if ( isView ) {
344  for (size_t i = 0; i < inputsHandle->size(); ++i) {
345  inputs_.push_back(inputsHandle->ptrAt(i));
346  }
347  } else {
348  iEvent.getByLabel(src_,pfinputsHandleAsFwdPtr);
349  for (size_t i = 0; i < pfinputsHandleAsFwdPtr->size(); ++i) {
350  if ( (*pfinputsHandleAsFwdPtr)[i].ptr().isAvailable() ) {
351  inputs_.push_back( (*pfinputsHandleAsFwdPtr)[i].ptr() );
352  }
353  else if ( (*pfinputsHandleAsFwdPtr)[i].backPtr().isAvailable() ) {
354  inputs_.push_back( (*pfinputsHandleAsFwdPtr)[i].backPtr() );
355  }
356  }
357  }
358  LogDebug("VirtualJetProducer") << "Got inputs\n";
359 
360  // Convert candidates to fastjet::PseudoJets.
361  // Also correct to Primary Vertex. Will modify fjInputs_
362  // and use inputs_
363  fjInputs_.reserve(inputs_.size());
364  inputTowers();
365  LogDebug("VirtualJetProducer") << "Inputted towers\n";
366 
367  // For Pileup subtraction using offset correction:
368  // Subtract pedestal.
369  if ( doPUOffsetCorr_ ) {
370  subtractor_->setDefinition(fjJetDefinition_);
372  subtractor_->calculatePedestal(fjInputs_);
373  subtractor_->subtractPedestal(fjInputs_);
374  LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
375  }
376  // Run algorithm. Will modify fjJets_ and allocate fjClusterSeq_.
377  // This will use fjInputs_
378  runAlgorithm( iEvent, iSetup );
379 
380  // if ( doPUOffsetCorr_ ) {
381  // subtractor_->setAlgorithm(fjClusterSeq_);
382  // }
383 
384  LogDebug("VirtualJetProducer") << "Ran algorithm\n";
385  // For Pileup subtraction using offset correction:
386  // Now we find jets and need to recalculate their energy,
387  // mark towers participated in jet,
388  // remove occupied towers from the list and recalculate mean and sigma
389  // put the initial towers collection to the jet,
390  // and subtract from initial towers in jet recalculated mean and sigma of towers
391  if ( doPUOffsetCorr_ ) {
392  LogDebug("VirtualJetProducer") << "Do PUOffsetCorr\n";
393  vector<fastjet::PseudoJet> orphanInput;
394  subtractor_->calculateOrphanInput(orphanInput);
395  subtractor_->calculatePedestal(orphanInput);
396  subtractor_->offsetCorrectJets();
397  }
398  // Write the output jets.
399  // This will (by default) call the member function template
400  // "writeJets", but can be overridden.
401  // this will use inputs_
402  output( iEvent, iSetup );
403  LogDebug("VirtualJetProducer") << "Wrote jets\n";
404 
405  return;
406 }
407 
408 //______________________________________________________________________________
409 
411 {
412  std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
413  inEnd = inputs_.end(), i = inBegin;
414  for (; i != inEnd; ++i ) {
416  if (edm::isNotFinite(input->pt())) continue;
417  if (input->et() <inputEtMin_) continue;
418  if (input->energy()<inputEMin_) continue;
419  if (isAnomalousTower(input)) continue;
420  if (input->pt() == 0) {
421  edm::LogError("NullTransverseMomentum") << "dropping input candidate with pt=0";
422  continue;
423  }
425  const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
427  fjInputs_.push_back(fastjet::PseudoJet(ct.px(),ct.py(),ct.pz(),ct.energy()));
428  //std::cout << "tower:" << *tower << '\n';
429  }
430  else {
431  /*
432  if(makePFJet(jetTypeE)) {
433  reco::PFCandidate* pfc = (reco::PFCandidate*)input.get();
434  std::cout << "PF cand:" << *pfc << '\n';
435  }
436  */
437  fjInputs_.push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
438  input->energy()));
439  }
440  fjInputs_.back().set_user_index(i - inBegin);
441  }
442 
443  if ( restrictInputs_ && fjInputs_.size() > maxInputs_ ) {
445  std::sort(fjInputs_.begin(), fjInputs_.end(), pTComparator);
446  fjInputs_.resize(maxInputs_);
447  edm::LogWarning("JetRecoTooManyEntries") << "Too many inputs in the event, limiting to first " << maxInputs_ << ". Output is suspect.";
448  }
449 }
450 
451 //______________________________________________________________________________
453 {
454  if (!makeCaloJet(jetTypeE))
455  return false;
456  else
457  return (*anomalousTowerDef_)(*input);
458 }
459 
460 //------------------------------------------------------------------------------
461 // This is pure virtual.
462 //______________________________________________________________________________
463 // void VirtualJetProducer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup,
464 // std::vector<edm::Ptr<reco::Candidate> > const & inputs_);
465 
466 //______________________________________________________________________________
467 void VirtualJetProducer::copyConstituents(const vector<fastjet::PseudoJet>& fjConstituents,
468  reco::Jet* jet)
469 {
470  for (unsigned int i=0;i<fjConstituents.size();++i) {
471  int index = fjConstituents[i].user_index();
472  if ( index >= 0 && static_cast<unsigned int>(index) < inputs_.size() )
473  jet->addDaughter(inputs_[index]);
474  }
475 }
476 
477 
478 //______________________________________________________________________________
479 vector<reco::CandidatePtr>
480 VirtualJetProducer::getConstituents(const vector<fastjet::PseudoJet>&fjConstituents)
481 {
482  vector<reco::CandidatePtr> result;
483  for (unsigned int i=0;i<fjConstituents.size();i++) {
484  int index = fjConstituents[i].user_index();
485  if ( index >= 0 && static_cast<unsigned int>(index) < inputs_.size() ) {
486  reco::CandidatePtr candidate = inputs_[index];
487  result.push_back(candidate);
488  }
489  }
490  return result;
491 }
492 
493 
494 //_____________________________________________________________________________
495 
497 {
498  // Write jets and constitutents. Will use fjJets_, inputs_
499  // and fjClusterSeq_
500 
501  if ( !writeCompound_ ) {
502  switch( jetTypeE ) {
503  case JetType::CaloJet :
504  writeJets<reco::CaloJet>( iEvent, iSetup);
505  break;
506  case JetType::PFJet :
507  writeJets<reco::PFJet>( iEvent, iSetup);
508  break;
509  case JetType::GenJet :
510  writeJets<reco::GenJet>( iEvent, iSetup);
511  break;
512  case JetType::TrackJet :
513  writeJets<reco::TrackJet>( iEvent, iSetup);
514  break;
515  case JetType::PFClusterJet :
516  writeJets<reco::PFClusterJet>( iEvent, iSetup);
517  break;
518  case JetType::BasicJet :
519  writeJets<reco::BasicJet>( iEvent, iSetup);
520  break;
521  default:
522  throw cms::Exception("InvalidInput") << "invalid jet type in VirtualJetProducer\n";
523  break;
524  };
525  } else {
526  // Write jets and constitutents.
527  switch( jetTypeE ) {
528  case JetType::CaloJet :
529  writeCompoundJets<reco::CaloJet>( iEvent, iSetup );
530  break;
531  case JetType::PFJet :
532  writeCompoundJets<reco::PFJet>( iEvent, iSetup );
533  break;
534  case JetType::GenJet :
535  writeCompoundJets<reco::GenJet>( iEvent, iSetup );
536  break;
537  case JetType::BasicJet :
538  writeCompoundJets<reco::BasicJet>( iEvent, iSetup );
539  break;
540  default:
541  throw cms::Exception("InvalidInput") << "invalid jet type in CompoundJetProducer\n";
542  break;
543  };
544  }
545 }
546 
547 namespace {
548 template< typename T >
549 struct Area { static float get(T const &) {return 0;}};
550 
551 template<>
552 struct Area<reco::CaloJet>{ static float get(reco::CaloJet const & jet) {
553  return jet.getSpecific().mTowersArea;
554 }
555 };
556 }
557 
558 template< typename T >
560 {
561  // std::cout << "writeJets " << typeid(T).name()
562  // << (doRhoFastjet_ ? " doRhoFastjet " : "")
563  // << (doAreaFastjet_ ? " doAreaFastjet " : "")
564  // << (doAreaDiskApprox_ ? " doAreaDiskApprox " : "")
565  // << std::endl;
566 
567  if (doRhoFastjet_) {
568  // declare jet collection without the two jets,
569  // for unbiased background estimation.
570  std::vector<fastjet::PseudoJet> fjexcluded_jets;
571  fjexcluded_jets=fjJets_;
572 
573  if(fjexcluded_jets.size()>2) fjexcluded_jets.resize(nExclude_);
574 
576  std::auto_ptr<std::vector<double> > rhos(new std::vector<double>);
577  std::auto_ptr<std::vector<double> > sigmas(new std::vector<double>);
578  int nEta = puCenters_.size();
579  rhos->reserve(nEta);
580  sigmas->reserve(nEta);
581  fastjet::ClusterSequenceAreaBase const* clusterSequenceWithArea =
582  dynamic_cast<fastjet::ClusterSequenceAreaBase const *> ( &*fjClusterSeq_ );
583 
584 
585  for(int ie = 0; ie < nEta; ++ie){
586  double eta = puCenters_[ie];
587  double etamin=eta-puWidth_;
588  double etamax=eta+puWidth_;
589  fastjet::RangeDefinition range_rho(etamin,etamax);
590  fastjet::BackgroundEstimator bkgestim(*clusterSequenceWithArea,range_rho);
591  bkgestim.set_excluded_jets(fjexcluded_jets);
592  rhos->push_back(bkgestim.rho());
593  sigmas->push_back(bkgestim.sigma());
594  }
595  iEvent.put(rhos,"rhos");
596  iEvent.put(sigmas,"sigmas");
597  }else{
598  std::auto_ptr<double> rho(new double(0.0));
599  std::auto_ptr<double> sigma(new double(0.0));
600  double mean_area = 0;
601 
602  fastjet::ClusterSequenceAreaBase const* clusterSequenceWithArea =
603  dynamic_cast<fastjet::ClusterSequenceAreaBase const *> ( &*fjClusterSeq_ );
604  /*
605  const double nemptyjets = clusterSequenceWithArea->n_empty_jets(*fjRangeDef_);
606  if(( nemptyjets < -15 ) || ( nemptyjets > fjRangeDef_->area()+ 15)) {
607  edm::LogWarning("StrangeNEmtpyJets") << "n_empty_jets is : " << clusterSequenceWithArea->n_empty_jets(*fjRangeDef_) << " with range " << fjRangeDef_->description() << ".";
608  }
609  */
610  clusterSequenceWithArea->get_median_rho_and_sigma(*fjRangeDef_,false,*rho,*sigma,mean_area);
611  if((*rho < 0)|| (edm::isNotFinite(*rho))) {
612  edm::LogError("BadRho") << "rho value is " << *rho << " area:" << mean_area << " and n_empty_jets: " << clusterSequenceWithArea->n_empty_jets(*fjRangeDef_) << " with range " << fjRangeDef_->description()
613  <<". Setting rho to rezo.";
614  *rho = 0;
615  }
616  iEvent.put(rho,"rho");
617  iEvent.put(sigma,"sigma");
618  }
619  } // doRhoFastjet_
620 
621  // produce output jet collection
622 
623  using namespace reco;
624 
625  std::auto_ptr<std::vector<T> > jets(new std::vector<T>() );
626  jets->reserve(fjJets_.size());
627 
628  // Distance between jet centers and overlap area -- for disk-based area calculation
629  using RIJ = std::pair<double,double>;
630  std::vector<std::vector<RIJ> > rij(fjJets_.size());
631 
632  float etaJ[fjJets_.size()], phiJ[fjJets_.size()];
633  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));};
634  for (auto ijet=0U;ijet<fjJets_.size();++ijet) {
635  float x = fjJets_[ijet].px();
636  float y = fjJets_[ijet].py();
637  float z = fjJets_[ijet].pz();
638  phiJ[ijet] = vdt::fast_atan2(y,x);
639  etaJ[ijet] =etaFromXYZ(x,y,z);
640  }
641  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
642  // allocate this jet
643  T jet;
644  // get the fastjet jet
645  const fastjet::PseudoJet& fjJet = fjJets_[ijet];
646  // get the constituents from fastjet
647  std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
648  // convert them to CandidatePtr vector
649  std::vector<CandidatePtr> constituents =
650  getConstituents(fjConstituents);
651 
652  // calcuate the jet area
653  double jetArea=0.0;
654  if ( doAreaFastjet_ && fjJet.has_area() ) {
655  jetArea = fjJet.area();
656  }
657  else if ( doAreaDiskApprox_ ) {
658  // Here it is assumed that fjJets_ is in decreasing order of pT,
659  // which should happen in FastjetJetProducer::runAlgorithm()
660  jetArea = M_PI;
661  if (0!=ijet) {
662  std::vector<RIJ>& distance = rij[ijet];
663  distance.resize(ijet);
664  for (unsigned jJet = 0; jJet < ijet; ++jJet) {
665  distance[jJet].first = std::sqrt(reco::deltaR2(etaJ[ijet],phiJ[ijet], etaJ[jJet],phiJ[jJet])) / rParam_;
666  distance[jJet].second = reco::helper::VirtualJetProducerHelper::intersection(distance[jJet].first);
667  jetArea -=distance[jJet].second;
668  for (unsigned kJet = 0; kJet < jJet; ++kJet) {
669  jetArea += reco::helper::VirtualJetProducerHelper::intersection(distance[jJet].first, distance[kJet].first, rij[jJet][kJet].first,
670  distance[jJet].second, distance[kJet].second, rij[jJet][kJet].second);
671  } // end loop over harder jets
672  } // end loop over harder jets
673  }
674  jetArea *= rParam_;
675  jetArea *= rParam_;
676  }
677  // write the specifics to the jet (simultaneously sets 4-vector, vertex).
678  // These are overridden functions that will call the appropriate
679  // specific allocator.
680  writeSpecific(jet,
681  Particle::LorentzVector(fjJet.px(),
682  fjJet.py(),
683  fjJet.pz(),
684  fjJet.E()),
685  vertex_,
686  constituents, iSetup);
687 
688  jet.setJetArea (jetArea);
689 
690  if(doPUOffsetCorr_){
691  jet.setPileup(subtractor_->getPileUpEnergy(ijet));
692  }else{
693  jet.setPileup (0.0);
694  }
695 
696 
697  // std::cout << "area " << ijet << " " << jetArea << " " << Area<T>::get(jet) << std::endl;
698 
699  // add to the list
700  jets->push_back(jet);
701  }
702  // put the jets in the collection
703  iEvent.put(jets,jetCollInstanceName_);
704 }
705 
707 template< class T>
709 {
710  if ( verbosity_ >= 1 ) {
711  std::cout << "<VirtualJetProducer::writeCompoundJets (moduleLabel = " << moduleLabel_ << ")>:" << std::endl;
712  }
713 
714  // get a list of output jets
715  std::auto_ptr<reco::BasicJetCollection> jetCollection( new reco::BasicJetCollection() );
716  // get a list of output subjets
717  std::auto_ptr<std::vector<T> > subjetCollection( new std::vector<T>() );
718 
719  // This will store the handle for the subjets after we write them
720  edm::OrphanHandle< std::vector<T> > subjetHandleAfterPut;
721  // this is the mapping of subjet to hard jet
722  std::vector< std::vector<int> > indices;
723  // this is the list of hardjet 4-momenta
724  std::vector<math::XYZTLorentzVector> p4_hardJets;
725  // this is the hardjet areas
726  std::vector<double> area_hardJets;
727 
728  // Loop over the hard jets
729  std::vector<fastjet::PseudoJet>::const_iterator it = fjJets_.begin(),
730  iEnd = fjJets_.end(),
731  iBegin = fjJets_.begin();
732  indices.resize( fjJets_.size() );
733  for ( ; it != iEnd; ++it ) {
734  fastjet::PseudoJet const & localJet = *it;
735  unsigned int jetIndex = it - iBegin;
736  // Get the 4-vector for the hard jet
737  p4_hardJets.push_back( math::XYZTLorentzVector(localJet.px(), localJet.py(), localJet.pz(), localJet.e() ));
738  double localJetArea = 0.0;
739  if ( doAreaFastjet_ && localJet.has_area() ) {
740  localJetArea = localJet.area();
741  }
742  area_hardJets.push_back( localJetArea );
743 
744  // create the subjet list
745  std::vector<fastjet::PseudoJet> constituents;
746  if ( it->has_pieces() ) {
747  constituents = it->pieces();
748  } else if ( it->has_constituents() ) {
749  constituents = it->constituents();
750  }
751 
752  std::vector<fastjet::PseudoJet>::const_iterator itSubJetBegin = constituents.begin(),
753  itSubJet = itSubJetBegin, itSubJetEnd = constituents.end();
754  for (; itSubJet != itSubJetEnd; ++itSubJet ){
755 
756  fastjet::PseudoJet const & subjet = *itSubJet;
757  if ( verbosity_ >= 1 ) {
758  std::cout << "subjet #" << (itSubJet - itSubJetBegin) << ": Pt = " << subjet.pt() << ", eta = " << subjet.eta() << ", phi = " << subjet.phi() << ", mass = " << subjet.m()
759  << " (#constituents = " << subjet.constituents().size() << ")" << std::endl;
760  std::vector<fastjet::PseudoJet> subjet_constituents = subjet.constituents();
761  int idx_constituent = 0;
762  for ( std::vector<fastjet::PseudoJet>::const_iterator constituent = subjet_constituents.begin();
763  constituent != subjet_constituents.end(); ++constituent ) {
764  if ( constituent->pt() < 1.e-3 ) continue; // CV: skip ghosts
765  std::cout << " constituent #" << idx_constituent << ": Pt = " << constituent->pt() << ", eta = " << constituent->eta() << ", phi = " << constituent->phi() << ","
766  << " mass = " << constituent->m() << std::endl;
767  ++idx_constituent;
768  }
769  }
770 
771  if ( verbosity_ >= 1 ) {
772  std::cout << "subjet #" << (itSubJet - itSubJetBegin) << ": Pt = " << subjet.pt() << ", eta = " << subjet.eta() << ", phi = " << subjet.phi() << ", mass = " << subjet.m()
773  << " (#constituents = " << subjet.constituents().size() << ")" << std::endl;
774  std::vector<fastjet::PseudoJet> subjet_constituents = subjet.constituents();
775  int idx_constituent = 0;
776  for ( std::vector<fastjet::PseudoJet>::const_iterator constituent = subjet_constituents.begin();
777  constituent != subjet_constituents.end(); ++constituent ) {
778  if ( constituent->pt() < 1.e-3 ) continue; // CV: skip ghosts
779  std::cout << " constituent #" << idx_constituent << ": Pt = " << constituent->pt() << ", eta = " << constituent->eta() << ", phi = " << constituent->phi() << ","
780  << " mass = " << constituent->m() << std::endl;
781  ++idx_constituent;
782  }
783  }
784 
785  math::XYZTLorentzVector p4Subjet(subjet.px(), subjet.py(), subjet.pz(), subjet.e() );
787 
788  // This will hold ptr's to the subjets
789  std::vector<reco::CandidatePtr> subjetConstituents;
790 
791  // Get the transient subjet constituents from fastjet
792  std::vector<fastjet::PseudoJet> subjetFastjetConstituents = subjet.constituents();
793  std::vector<reco::CandidatePtr> constituents =
794  getConstituents(subjetFastjetConstituents );
795 
796  indices[jetIndex].push_back( subjetCollection->size() );
797 
798  // Add the concrete subjet type to the subjet list to write to event record
799  T jet;
800  reco::writeSpecific( jet, p4Subjet, point, constituents, iSetup);
801  double subjetArea = 0.0;
802  if ( doAreaFastjet_ && itSubJet->has_area() ){
803  subjetArea = itSubJet->area();
804  }
805  jet.setJetArea( subjetArea );
806  subjetCollection->push_back( jet );
807  }
808  }
809 
810  // put subjets into event record
811  subjetHandleAfterPut = iEvent.put( subjetCollection, jetCollInstanceName_ );
812 
813  // Now create the hard jets with ptr's to the subjets as constituents
814  std::vector<math::XYZTLorentzVector>::const_iterator ip4 = p4_hardJets.begin(),
815  ip4Begin = p4_hardJets.begin(),
816  ip4End = p4_hardJets.end();
817 
818  for ( ; ip4 != ip4End; ++ip4 ) {
819  int p4_index = ip4 - ip4Begin;
820  std::vector<int> & ind = indices[p4_index];
821  std::vector<reco::CandidatePtr> i_hardJetConstituents;
822  // Add the subjets to the hard jet
823  for( std::vector<int>::const_iterator isub = ind.begin();
824  isub != ind.end(); ++isub ) {
825  reco::CandidatePtr candPtr( subjetHandleAfterPut, *isub, false );
826  i_hardJetConstituents.push_back( candPtr );
827  }
829  reco::BasicJet toput( *ip4, point, i_hardJetConstituents);
830  toput.setJetArea( area_hardJets[ip4 - ip4Begin] );
831  jetCollection->push_back( toput );
832  }
833 
834  // put hard jets into event record
835  iEvent.put( jetCollection);
836 }
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:42
boost::shared_ptr< fastjet::AreaDefinition > AreaDefinitionPtr
void writeJets(edm::Event &iEvent, edm::EventSetup const &iSetup)
JetType::Type jetTypeE
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
T getUntrackedParameter(std::string const &, T const &) const
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:105
virtual std::vector< reco::CandidatePtr > getConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents)
std::vector< fastjet::PseudoJet > fjJets_
virtual void inputTowers()
Base class for all types of Jets.
Definition: Jet.h:20
boost::shared_ptr< fastjet::JetDefinition::Plugin > PluginPtr
Definition: DDAxes.h:10
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::string puSubtractorName_
std::vector< double > puCenters_
T eta() const
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
float float float z
virtual void runAlgorithm(edm::Event &iEvent, const edm::EventSetup &iSetup)=0
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)
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
static const char * names[]
std::vector< fastjet::PseudoJet > fjInputs_
int iEvent
Definition: GenABIO.cc:243
std::string jetCollInstanceName_
boost::shared_ptr< PileUpSubtractor > subtractor_
bool isNotFinite(T x)
Definition: isFinite.h:10
std::vector< edm::Ptr< reco::Candidate > > inputs_
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:143
const T & max(const T &a, const T &b)
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:18
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
ClusterSequencePtr fjClusterSeq_
virtual void makeProduces(std::string s, std::string tag="")
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:31
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_
bool first
Definition: L1TdeRCT.cc:79
ActiveAreaSpecPtr fjActiveArea_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
#define M_PI
Definition: BFit3D.cc:3
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
double rho()
synonym for median rho [[do we have this? Both?]]
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:58
double sigma()
get the sigma
edm::EventID id() const
Definition: EventBase.h:56
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
Definition: DDAxes.h:10
moduleLabel_(iConfig.getParameter< string >("@module_label"))
long double T
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
bool makeCaloJet(const JetType::Type &fTag)
static Type byName(const std::string &name)
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::XYZPoint Point
point in the space
Definition: LeafCandidate.h:30
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