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