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