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 
20 
34 
35 //#include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
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 
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 
82 void VirtualJetProducer::makeProduces( std::string alias, std::string tag )
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 {
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");
190  else puSubtractorName_ = std::string();
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 
276  produces<std::vector<double> >("rhos");
277  produces<std::vector<double> >("sigmas");
278  produces<double>("rho");
279  produces<double>("sigma");
280 
281 
282 }
283 
284 //______________________________________________________________________________
286 {
287 }
288 
289 
291 // implementation of member functions
293 
294 //______________________________________________________________________________
296 {
297 
298  // If requested, set the fastjet random seed to a deterministic function
299  // of the run/lumi/event.
300  // NOTE!!! The fastjet random number sequence is a global singleton.
301  // Thus, we have to create an object and get access to the global singleton
302  // in order to change it.
303  if ( useDeterministicSeed_ ) {
304  fastjet::GhostedAreaSpec gas;
305  std::vector<int> seeds(2);
306  seeds[0] = std::max(iEvent.id().run(),minSeed_ + 3) + 3 * iEvent.id().event();
307  seeds[1] = std::max(iEvent.id().run(),minSeed_ + 5) + 5 * iEvent.id().event();
308  gas.set_random_status(seeds);
309  }
310 
311  LogDebug("VirtualJetProducer") << "Entered produce\n";
312  //determine signal vertex2
313  vertex_=reco::Jet::Point(0,0,0);
315  LogDebug("VirtualJetProducer") << "Adding PV info\n";
317  iEvent.getByLabel(srcPVs_,pvCollection);
318  if (pvCollection->size()>0) vertex_=pvCollection->begin()->position();
319  }
320 
321  // For Pileup subtraction using offset correction:
322  // set up geometry map
323  if ( doPUOffsetCorr_ ) {
324  subtractor_->setupGeometryMap(iEvent, iSetup);
325  }
326 
327  // clear data
328  LogDebug("VirtualJetProducer") << "Clear data\n";
329  fjInputs_.clear();
330  fjJets_.clear();
331  inputs_.clear();
332 
333  // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
335  iEvent.getByLabel(src_,inputsHandle);
336  for (size_t i = 0; i < inputsHandle->size(); ++i) {
337  inputs_.push_back(inputsHandle->ptrAt(i));
338  }
339  LogDebug("VirtualJetProducer") << "Got inputs\n";
340 
341  // Convert candidates to fastjet::PseudoJets.
342  // Also correct to Primary Vertex. Will modify fjInputs_
343  // and use inputs_
344  fjInputs_.reserve(inputs_.size());
345  inputTowers();
346  LogDebug("VirtualJetProducer") << "Inputted towers\n";
347 
348  // For Pileup subtraction using offset correction:
349  // Subtract pedestal.
350  if ( doPUOffsetCorr_ ) {
351  subtractor_->setDefinition(fjJetDefinition_);
353  subtractor_->calculatePedestal(fjInputs_);
354  subtractor_->subtractPedestal(fjInputs_);
355  LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
356  }
357  // Run algorithm. Will modify fjJets_ and allocate fjClusterSeq_.
358  // This will use fjInputs_
359  runAlgorithm( iEvent, iSetup );
360 
361  // if ( doPUOffsetCorr_ ) {
362  // subtractor_->setAlgorithm(fjClusterSeq_);
363  // }
364 
365  LogDebug("VirtualJetProducer") << "Ran algorithm\n";
366  // For Pileup subtraction using offset correction:
367  // Now we find jets and need to recalculate their energy,
368  // mark towers participated in jet,
369  // remove occupied towers from the list and recalculate mean and sigma
370  // put the initial towers collection to the jet,
371  // and subtract from initial towers in jet recalculated mean and sigma of towers
372  if ( doPUOffsetCorr_ ) {
373  LogDebug("VirtualJetProducer") << "Do PUOffsetCorr\n";
374  vector<fastjet::PseudoJet> orphanInput;
375  subtractor_->calculateOrphanInput(orphanInput);
376  subtractor_->calculatePedestal(orphanInput);
377  subtractor_->offsetCorrectJets();
378  }
379  // Write the output jets.
380  // This will (by default) call the member function template
381  // "writeJets", but can be overridden.
382  // this will use inputs_
383  output( iEvent, iSetup );
384  LogDebug("VirtualJetProducer") << "Wrote jets\n";
385 
386  return;
387 }
388 
389 //______________________________________________________________________________
390 
392 {
393  std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
394  inEnd = inputs_.end(), i = inBegin;
395  for (; i != inEnd; ++i ) {
397  if (std::isnan(input->pt())) continue;
398  if (input->et() <inputEtMin_) continue;
399  if (input->energy()<inputEMin_) continue;
400  if (isAnomalousTower(input)) continue;
401  if (input->pt() == 0) {
402  edm::LogError("NullTransverseMomentum") << "dropping input candidate with pt=0";
403  continue;
404  }
406  const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
408  fjInputs_.push_back(fastjet::PseudoJet(ct.px(),ct.py(),ct.pz(),ct.energy()));
409  //std::cout << "tower:" << *tower << '\n';
410  }
411  else {
412  /*
413  if(makePFJet(jetTypeE)) {
414  reco::PFCandidate* pfc = (reco::PFCandidate*)input.get();
415  std::cout << "PF cand:" << *pfc << '\n';
416  }
417  */
418  fjInputs_.push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
419  input->energy()));
420  }
421  fjInputs_.back().set_user_index(i - inBegin);
422  }
423 
424  if ( restrictInputs_ && fjInputs_.size() > maxInputs_ ) {
426  std::sort(fjInputs_.begin(), fjInputs_.end(), pTComparator);
427  fjInputs_.resize(maxInputs_);
428  edm::LogWarning("JetRecoTooManyEntries") << "Too many inputs in the event, limiting to first " << maxInputs_ << ". Output is suspect.";
429  }
430 }
431 
432 //______________________________________________________________________________
434 {
435  if (!makeCaloJet(jetTypeE))
436  return false;
437  else
438  return (*anomalousTowerDef_)(*input);
439 }
440 
441 //------------------------------------------------------------------------------
442 // This is pure virtual.
443 //______________________________________________________________________________
444 // void VirtualJetProducer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup,
445 // std::vector<edm::Ptr<reco::Candidate> > const & inputs_);
446 
447 //______________________________________________________________________________
448 void VirtualJetProducer::copyConstituents(const vector<fastjet::PseudoJet>& fjConstituents,
449  reco::Jet* jet)
450 {
451  for (unsigned int i=0;i<fjConstituents.size();++i) {
452  int index = fjConstituents[i].user_index();
453  if ( index >= 0 && static_cast<unsigned int>(index) < inputs_.size() )
454  jet->addDaughter(inputs_[index]);
455  }
456 }
457 
458 
459 //______________________________________________________________________________
460 vector<reco::CandidatePtr>
461 VirtualJetProducer::getConstituents(const vector<fastjet::PseudoJet>&fjConstituents)
462 {
463  vector<reco::CandidatePtr> result;
464  for (unsigned int i=0;i<fjConstituents.size();i++) {
465  int index = fjConstituents[i].user_index();
466  if ( index >= 0 && static_cast<unsigned int>(index) < inputs_.size() ) {
467  reco::CandidatePtr candidate = inputs_[index];
468  result.push_back(candidate);
469  }
470  }
471  return result;
472 }
473 
474 
475 //_____________________________________________________________________________
476 
478 {
479  // Write jets and constitutents. Will use fjJets_, inputs_
480  // and fjClusterSeq_
481 
482  if ( !writeCompound_ ) {
483  switch( jetTypeE ) {
484  case JetType::CaloJet :
485  writeJets<reco::CaloJet>( iEvent, iSetup);
486  break;
487  case JetType::PFJet :
488  writeJets<reco::PFJet>( iEvent, iSetup);
489  break;
490  case JetType::GenJet :
491  writeJets<reco::GenJet>( iEvent, iSetup);
492  break;
493  case JetType::TrackJet :
494  writeJets<reco::TrackJet>( iEvent, iSetup);
495  break;
496  case JetType::PFClusterJet :
497  writeJets<reco::PFClusterJet>( iEvent, iSetup);
498  break;
499  case JetType::BasicJet :
500  writeJets<reco::BasicJet>( iEvent, iSetup);
501  break;
502  default:
503  throw cms::Exception("InvalidInput") << "invalid jet type in VirtualJetProducer\n";
504  break;
505  };
506  } else {
507  // Write jets and constitutents.
508  switch( jetTypeE ) {
509  case JetType::CaloJet :
510  writeCompoundJets<reco::CaloJet>( iEvent, iSetup );
511  break;
512  case JetType::PFJet :
513  writeCompoundJets<reco::PFJet>( iEvent, iSetup );
514  break;
515  case JetType::GenJet :
516  writeCompoundJets<reco::GenJet>( iEvent, iSetup );
517  break;
518  case JetType::BasicJet :
519  writeCompoundJets<reco::BasicJet>( iEvent, iSetup );
520  break;
521  default:
522  throw cms::Exception("InvalidInput") << "invalid jet type in CompoundJetProducer\n";
523  break;
524  };
525  }
526 
527 }
528 
529 template< typename T >
531 {
532  if (doRhoFastjet_) {
533  // declare jet collection without the two jets,
534  // for unbiased background estimation.
535  std::vector<fastjet::PseudoJet> fjexcluded_jets;
536  fjexcluded_jets=fjJets_;
537 
538  if(fjexcluded_jets.size()>2) fjexcluded_jets.resize(nExclude_);
539 
541  std::auto_ptr<std::vector<double> > rhos(new std::vector<double>);
542  std::auto_ptr<std::vector<double> > sigmas(new std::vector<double>);
543  int nEta = puCenters_.size();
544  rhos->reserve(nEta);
545  sigmas->reserve(nEta);
546  fastjet::ClusterSequenceAreaBase const* clusterSequenceWithArea =
547  dynamic_cast<fastjet::ClusterSequenceAreaBase const *> ( &*fjClusterSeq_ );
548 
549 
550  for(int ie = 0; ie < nEta; ++ie){
551  double eta = puCenters_[ie];
552  double etamin=eta-puWidth_;
553  double etamax=eta+puWidth_;
554  fastjet::RangeDefinition range_rho(etamin,etamax);
555  fastjet::BackgroundEstimator bkgestim(*clusterSequenceWithArea,range_rho);
556  bkgestim.set_excluded_jets(fjexcluded_jets);
557  rhos->push_back(bkgestim.rho());
558  sigmas->push_back(bkgestim.sigma());
559  }
560  iEvent.put(rhos,"rhos");
561  iEvent.put(sigmas,"sigmas");
562  }else{
563  std::auto_ptr<double> rho(new double(0.0));
564  std::auto_ptr<double> sigma(new double(0.0));
565  double mean_area = 0;
566 
567  fastjet::ClusterSequenceAreaBase const* clusterSequenceWithArea =
568  dynamic_cast<fastjet::ClusterSequenceAreaBase const *> ( &*fjClusterSeq_ );
569  /*
570  const double nemptyjets = clusterSequenceWithArea->n_empty_jets(*fjRangeDef_);
571  if(( nemptyjets < -15 ) || ( nemptyjets > fjRangeDef_->area()+ 15)) {
572  edm::LogWarning("StrangeNEmtpyJets") << "n_empty_jets is : " << clusterSequenceWithArea->n_empty_jets(*fjRangeDef_) << " with range " << fjRangeDef_->description() << ".";
573  }
574  */
575  clusterSequenceWithArea->get_median_rho_and_sigma(*fjRangeDef_,false,*rho,*sigma,mean_area);
576  if((*rho < 0)|| (std::isnan(*rho))) {
577  edm::LogError("BadRho") << "rho value is " << *rho << " area:" << mean_area << " and n_empty_jets: " << clusterSequenceWithArea->n_empty_jets(*fjRangeDef_) << " with range " << fjRangeDef_->description()
578  <<". Setting rho to rezo.";
579  *rho = 0;
580  }
581  iEvent.put(rho,"rho");
582  iEvent.put(sigma,"sigma");
583  }
584  } // doRhoFastjet_
585 
586  // produce output jet collection
587 
588  using namespace reco;
589 
590  std::auto_ptr<std::vector<T> > jets(new std::vector<T>() );
591  jets->reserve(fjJets_.size());
592 
593  // Distance between jet centers -- for disk-based area calculation
594  std::vector<std::vector<double> > rij(fjJets_.size());
595 
596  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
597  // allocate this jet
598  T jet;
599  // get the fastjet jet
600  const fastjet::PseudoJet& fjJet = fjJets_[ijet];
601  // get the constituents from fastjet
602  std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
603  // convert them to CandidatePtr vector
604  std::vector<CandidatePtr> constituents =
605  getConstituents(fjConstituents);
606 
607 
608  // calcuate the jet area
609  double jetArea=0.0;
610  if ( doAreaFastjet_ && fjJet.has_area() ) {
611  jetArea = fjJet.area();
612  }
613  else if ( doAreaDiskApprox_ ) {
614  // Here it is assumed that fjJets_ is in decreasing order of pT,
615  // which should happen in FastjetJetProducer::runAlgorithm()
616  jetArea = M_PI;
617  if (ijet) {
618  std::vector<double>& distance = rij[ijet];
619  distance.resize(ijet);
620  for (unsigned jJet = 0; jJet < ijet; ++jJet) {
621  distance[jJet] = reco::deltaR(fjJets_[ijet],fjJets_[jJet]) / rParam_;
622  jetArea -= reco::helper::VirtualJetProducerHelper::intersection(distance[jJet]);
623  for (unsigned kJet = 0; kJet < jJet; ++kJet) {
624  jetArea += reco::helper::VirtualJetProducerHelper::intersection(distance[jJet], distance[kJet], rij[jJet][kJet]);
625  } // end loop over harder jets
626  } // end loop over harder jets
627  }
628  jetArea *= rParam_;
629  jetArea *= rParam_;
630  }
631  // write the specifics to the jet (simultaneously sets 4-vector, vertex).
632  // These are overridden functions that will call the appropriate
633  // specific allocator.
634  writeSpecific(jet,
635  Particle::LorentzVector(fjJet.px(),
636  fjJet.py(),
637  fjJet.pz(),
638  fjJet.E()),
639  vertex_,
640  constituents, iSetup);
641 
642  jet.setJetArea (jetArea);
643 
644  if(doPUOffsetCorr_){
645  jet.setPileup(subtractor_->getPileUpEnergy(ijet));
646  }else{
647  jet.setPileup (0.0);
648  }
649 
650  // add to the list
651  jets->push_back(jet);
652  }
653  // put the jets in the collection
654  iEvent.put(jets,jetCollInstanceName_);
655 
656 
657 }
658 
659 
660 
662 template< class T>
664 {
665 
666  // get a list of output jets
667  std::auto_ptr<reco::BasicJetCollection> jetCollection( new reco::BasicJetCollection() );
668  // get a list of output subjets
669  std::auto_ptr<std::vector<T> > subjetCollection( new std::vector<T>() );
670 
671  // This will store the handle for the subjets after we write them
672  edm::OrphanHandle< std::vector<T> > subjetHandleAfterPut;
673  // this is the mapping of subjet to hard jet
674  std::vector< std::vector<int> > indices;
675  // this is the list of hardjet 4-momenta
676  std::vector<math::XYZTLorentzVector> p4_hardJets;
677  // this is the hardjet areas
678  std::vector<double> area_hardJets;
679 
680 
681  // Loop over the hard jets
682  std::vector<fastjet::PseudoJet>::const_iterator it = fjJets_.begin(),
683  iEnd = fjJets_.end(),
684  iBegin = fjJets_.begin();
685  indices.resize( fjJets_.size() );
686  for ( ; it != iEnd; ++it ) {
687  fastjet::PseudoJet const & localJet = *it;
688  unsigned int jetIndex = it - iBegin;
689  // Get the 4-vector for the hard jet
690  p4_hardJets.push_back( math::XYZTLorentzVector(localJet.px(), localJet.py(), localJet.pz(), localJet.e() ));
691  double localJetArea = 0.0;
692  if ( doAreaFastjet_ && localJet.has_area() ) {
693  localJetArea = localJet.area();
694  }
695  area_hardJets.push_back( localJetArea );
696 
697  // create the subjet list
698  std::vector<fastjet::PseudoJet> constituents;
699  if ( it->has_pieces() ) {
700  constituents = it->pieces();
701  } else {
702  constituents=it->constituents();
703  }
704 
705 
706  std::vector<fastjet::PseudoJet>::const_iterator itSubJetBegin = constituents.begin(),
707  itSubJet = itSubJetBegin, itSubJetEnd = constituents.end();
708  for (; itSubJet != itSubJetEnd; ++itSubJet ){
709 
710  fastjet::PseudoJet const & subjet = *itSubJet;
711 
712  math::XYZTLorentzVector p4Subjet(subjet.px(), subjet.py(), subjet.pz(), subjet.e() );
714 
715  // This will hold ptr's to the subjets
716  std::vector<reco::CandidatePtr> subjetConstituents;
717 
718  // Get the transient subjet constituents from fastjet
719  std::vector<fastjet::PseudoJet> subjetFastjetConstituents = subjet.constituents();
720  std::vector<reco::CandidatePtr> constituents =
721  getConstituents(subjetFastjetConstituents );
722 
723  indices[jetIndex].push_back( subjetCollection->size() );
724 
725  // Add the concrete subjet type to the subjet list to write to event record
726  T jet;
727  reco::writeSpecific( jet, p4Subjet, point, constituents, iSetup);
728  double subjetArea = 0.0;
729  if ( doAreaFastjet_ && itSubJet->has_area() ){
730  subjetArea = itSubJet->area();
731  }
732  jet.setJetArea( subjetArea );
733  subjetCollection->push_back( jet );
734 
735  }
736  }
737  // put subjets into event record
738  subjetHandleAfterPut = iEvent.put( subjetCollection, jetCollInstanceName_ );
739 
740 
741  // Now create the hard jets with ptr's to the subjets as constituents
742  std::vector<math::XYZTLorentzVector>::const_iterator ip4 = p4_hardJets.begin(),
743  ip4Begin = p4_hardJets.begin(),
744  ip4End = p4_hardJets.end();
745 
746  for ( ; ip4 != ip4End; ++ip4 ) {
747  int p4_index = ip4 - ip4Begin;
748  std::vector<int> & ind = indices[p4_index];
749  std::vector<reco::CandidatePtr> i_hardJetConstituents;
750  // Add the subjets to the hard jet
751  for( std::vector<int>::const_iterator isub = ind.begin();
752  isub != ind.end(); ++isub ) {
753  reco::CandidatePtr candPtr( subjetHandleAfterPut, *isub, false );
754  i_hardJetConstituents.push_back( candPtr );
755  }
757  reco::BasicJet toput( *ip4, point, i_hardJetConstituents);
758  toput.setJetArea( area_hardJets[ip4 - ip4Begin] );
759  jetCollection->push_back( toput );
760  }
761 
762  // put hard jets into event record
763  iEvent.put( jetCollection);
764 
765 }
#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_
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:21
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:21
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
double deltaR(double eta1, double phi1, double eta2, double phi2)
Definition: deltaR.h:19
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
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:26
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
virtual bool isAnomalousTower(reco::CandidatePtr input)
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:104
static const char * names[]
std::vector< fastjet::PseudoJet > fjInputs_
int iEvent
Definition: GenABIO.cc:243
std::string jetCollInstanceName_
boost::shared_ptr< PileUpSubtractor > subtractor_
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)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
bool isnan(float x)
Definition: math.h:13
ClusterSequencePtr fjClusterSeq_
virtual void makeProduces(std::string s, std::string tag="")
vector< PseudoJet > jets
tuple result
Definition: query.py:137
math::XYZPoint Point
point in the space
Definition: Particle.h:29
void set_excluded_jets(const std::vector< PseudoJet > &excluded_jets)
boost::shared_ptr< fastjet::JetDefinition > JetDefPtr
std::auto_ptr< AnomalousTower > anomalousTowerDef_
ActiveAreaSpecPtr fjActiveArea_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
#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?]]
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
AreaDefinitionPtr fjAreaDefinition_
virtual void output(edm::Event &iEvent, edm::EventSetup const &iSetup)
long double T
static const HistoName names[]
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25
bool makeCaloJet(const JetType::Type &fTag)
static Type byName(const std::string &name)
T get(const Candidate &c)
Definition: component.h:56
*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:29
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