CMS 3D CMS Logo

Classes | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes

VirtualJetProducer Class Reference

#include <VirtualJetProducer.h>

Inheritance diagram for VirtualJetProducer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper cms::CompoundJetProducer cms::SubEventGenJetProducer CMSInsideOutJetProducer FastjetJetProducer SubjetFilterJetProducer cms::CATopJetProducer cms::SubJetProducer

List of all members.

Classes

struct  JetType

Public Types

typedef boost::shared_ptr
< fastjet::ActiveAreaSpec > 
ActiveAreaSpecPtr
typedef boost::shared_ptr
< fastjet::ClusterSequence > 
ClusterSequencePtr
typedef boost::shared_ptr
< fastjet::JetDefinition > 
JetDefPtr
typedef boost::shared_ptr
< fastjet::JetDefinition::Plugin > 
PluginPtr
typedef boost::shared_ptr
< fastjet::RangeDefinition > 
RangeDefPtr

Public Member Functions

std::string jetType () const
virtual void produce (edm::Event &iEvent, const edm::EventSetup &iSetup)
 VirtualJetProducer (const edm::ParameterSet &iConfig)
virtual ~VirtualJetProducer ()

Protected Member Functions

virtual void copyConstituents (const std::vector< fastjet::PseudoJet > &fjConstituents, reco::Jet *jet)
virtual std::vector
< reco::CandidatePtr
getConstituents (const std::vector< fastjet::PseudoJet > &fjConstituents)
virtual void inputTowers ()
virtual bool isAnomalousTower (reco::CandidatePtr input)
bool makeBasicJet (const JetType::Type &fTag)
bool makeCaloJet (const JetType::Type &fTag)
bool makeGenJet (const JetType::Type &fTag)
bool makePFJet (const JetType::Type &fTag)
virtual void makeProduces (std::string s, std::string tag="")
bool makeTrackJet (const JetType::Type &fTag)
void offsetCorrectJets (std::vector< fastjet::PseudoJet > &orphanInput)
virtual void output (edm::Event &iEvent, edm::EventSetup const &iSetup)
virtual void runAlgorithm (edm::Event &iEvent, const edm::EventSetup &iSetup)=0
template<typename T >
void writeJets (edm::Event &iEvent, edm::EventSetup const &iSetup)

Protected Attributes

bool doAreaFastjet_
bool doFastJetNonUniform_
bool doPUOffsetCorr_
bool doPVCorrection_
bool doRhoFastjet_
ActiveAreaSpecPtr fjActiveArea_
ClusterSequencePtr fjClusterSeq_
std::vector< fastjet::PseudoJet > fjInputs_
JetDefPtr fjJetDefinition_
std::vector< fastjet::PseudoJet > fjJets_
PluginPtr fjPlugin_
RangeDefPtr fjRangeDef_
double inputEMin_
double inputEtMin_
std::vector< edm::Ptr
< reco::Candidate > > 
inputs_
std::string jetAlgorithm_
std::string jetCollInstanceName_
double jetPtMin_
std::string jetType_
JetType::Type jetTypeE
unsigned int maxBadEcalCells_
unsigned int maxBadHcalCells_
unsigned int maxInputs_
unsigned int maxProblematicEcalCells_
unsigned int maxProblematicHcalCells_
unsigned int maxRecoveredEcalCells_
unsigned int maxRecoveredHcalCells_
std::string moduleLabel_
std::vector< double > puCenters_
std::string puSubtractorName_
double puWidth_
bool restrictInputs_
double rParam_
edm::InputTag src_
edm::InputTag srcPVs_
boost::shared_ptr
< PileUpSubtractor
subtractor_
reco::Particle::Point vertex_

Detailed Description

Definition at line 27 of file VirtualJetProducer.h.


Member Typedef Documentation

typedef boost::shared_ptr<fastjet::ActiveAreaSpec> VirtualJetProducer::ActiveAreaSpecPtr

Definition at line 76 of file VirtualJetProducer.h.

typedef boost::shared_ptr<fastjet::ClusterSequence> VirtualJetProducer::ClusterSequencePtr

Definition at line 73 of file VirtualJetProducer.h.

typedef boost::shared_ptr<fastjet::JetDefinition> VirtualJetProducer::JetDefPtr

Definition at line 75 of file VirtualJetProducer.h.

typedef boost::shared_ptr<fastjet::JetDefinition::Plugin> VirtualJetProducer::PluginPtr

Definition at line 74 of file VirtualJetProducer.h.

typedef boost::shared_ptr<fastjet::RangeDefinition> VirtualJetProducer::RangeDefPtr

Definition at line 77 of file VirtualJetProducer.h.


Constructor & Destructor Documentation

VirtualJetProducer::VirtualJetProducer ( const edm::ParameterSet iConfig) [explicit]

Definition at line 99 of file VirtualJetProducer.cc.

References ExpressReco_HICollisions_FallBack::alias, VirtualJetProducer::JetType::byName(), VirtualJetProducer::JetType::CaloJet, doAreaFastjet_, doFastJetNonUniform_, doPUOffsetCorr_, doRhoFastjet_, Exception, edm::ParameterSet::exists(), fjActiveArea_, fjJetDefinition_, fjPlugin_, fjRangeDef_, reco::get(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), jetAlgorithm_, jetCollInstanceName_, jetType_, jetTypeE, makeProduces(), maxInputs_, moduleLabel_, puCenters_, puSubtractorName_, puWidth_, restrictInputs_, rParam_, and subtractor_.

  : moduleLabel_   (iConfig.getParameter<string>       ("@module_label"))
  , src_           (iConfig.getParameter<edm::InputTag>("src"))
  , srcPVs_        (iConfig.getParameter<edm::InputTag>("srcPVs"))
  , jetType_       (iConfig.getParameter<string>       ("jetType"))
  , jetAlgorithm_  (iConfig.getParameter<string>       ("jetAlgorithm"))
  , rParam_        (iConfig.getParameter<double>       ("rParam"))
  , inputEtMin_    (iConfig.getParameter<double>       ("inputEtMin"))
  , inputEMin_     (iConfig.getParameter<double>       ("inputEMin"))
  , jetPtMin_      (iConfig.getParameter<double>       ("jetPtMin"))
  , doPVCorrection_(iConfig.getParameter<bool>         ("doPVCorrection"))
  , restrictInputs_(false)
  , maxInputs_(99999999)
  , doAreaFastjet_ (iConfig.getParameter<bool>         ("doAreaFastjet"))
  , doRhoFastjet_  (iConfig.getParameter<bool>         ("doRhoFastjet"))
  , doPUOffsetCorr_(iConfig.getParameter<bool>         ("doPUOffsetCorr"))
  , maxBadEcalCells_        (iConfig.getParameter<unsigned>("maxBadEcalCells"))
  , maxRecoveredEcalCells_  (iConfig.getParameter<unsigned>("maxRecoveredEcalCells"))
  , maxProblematicEcalCells_(iConfig.getParameter<unsigned>("maxProblematicEcalCells"))
  , maxBadHcalCells_        (iConfig.getParameter<unsigned>("maxBadHcalCells"))
  , maxRecoveredHcalCells_  (iConfig.getParameter<unsigned>("maxRecoveredHcalCells"))
  , maxProblematicHcalCells_(iConfig.getParameter<unsigned>("maxProblematicHcalCells"))
  , jetCollInstanceName_ ("")
{

  //
  // additional parameters to think about:
  // - overlap threshold (set to 0.75 for the time being)
  // - p parameter for generalized kT (set to -2 for the time being)
  // - fastjet PU subtraction parameters (not yet considered)
  //
  if (jetAlgorithm_=="SISCone") {
    fjPlugin_ = PluginPtr( new fastjet::SISConePlugin(rParam_,0.75,0,0.0,false,
                                                      fastjet::SISConePlugin::SM_pttilde) );
    fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(&*fjPlugin_) );
  }
  else if (jetAlgorithm_=="IterativeCone") {
    fjPlugin_ = PluginPtr(new fastjet::CMSIterativeConePlugin(rParam_,1.0));
    fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
  }
  else if (jetAlgorithm_=="CDFMidPoint") {
    fjPlugin_ = PluginPtr(new fastjet::CDFMidPointPlugin(rParam_,0.75));
    fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
  }
  else if (jetAlgorithm_=="ATLASCone") {
    fjPlugin_ = PluginPtr(new fastjet::ATLASConePlugin(rParam_));
    fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
  }
  else if (jetAlgorithm_=="Kt")
    fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::kt_algorithm,rParam_));
  else if (jetAlgorithm_=="CambridgeAachen")
    fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::cambridge_algorithm,
                                                           rParam_) );
  else if (jetAlgorithm_=="AntiKt")
    fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::antikt_algorithm,rParam_) );
  else if (jetAlgorithm_=="GeneralizedKt")
    fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::genkt_algorithm,
                                                            rParam_,-2) );
  else
    throw cms::Exception("Invalid jetAlgorithm")
      <<"Jet algorithm for VirtualJetProducer is invalid, Abort!\n";
  
  jetTypeE=JetType::byName(jetType_);

  if ( iConfig.exists("jetCollInstanceName") ) {
    jetCollInstanceName_ = iConfig.getParameter<string>("jetCollInstanceName");
  }

  if ( doPUOffsetCorr_ ) {
     if ( jetTypeE != JetType::CaloJet ) {
        throw cms::Exception("InvalidInput") << "Can only offset correct jets of type CaloJet";
     }
     
     if(iConfig.exists("subtractorName")) puSubtractorName_  =  iConfig.getParameter<string> ("subtractorName");
     else puSubtractorName_ = std::string();
     
     if(puSubtractorName_.empty()){
       edm::LogWarning("VirtualJetProducer") << "Pile Up correction on; however, pile up type is not specified. Using default... \n";
       subtractor_ =  boost::shared_ptr<PileUpSubtractor>(new PileUpSubtractor(iConfig));
     }else{
       subtractor_ =  boost::shared_ptr<PileUpSubtractor>(PileUpSubtractorFactory::get()->create( puSubtractorName_, iConfig));
     }
  }

  // do fasjet area / rho calcluation? => accept corresponding parameters
  if ( doAreaFastjet_ || doRhoFastjet_ ) {
    // Eta range of jets to be considered for Rho calculation
    // Should be at most (jet acceptance - jet radius)
    double rhoEtaMax=iConfig.getParameter<double>("Rho_EtaMax");
    // default Ghost_EtaMax should be 5
    double ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
    // default Active_Area_Repeats 1
    int    activeAreaRepeats = iConfig.getParameter<int> ("Active_Area_Repeats");
    // default GhostArea 0.01
    double ghostArea = iConfig.getParameter<double> ("GhostArea");
    fjActiveArea_ =  ActiveAreaSpecPtr(new fastjet::ActiveAreaSpec(ghostEtaMax,
                                                                   activeAreaRepeats,
                                                                   ghostArea));
    fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(rhoEtaMax) );
  } 

  // restrict inputs to first "maxInputs" towers?
  if ( iConfig.exists("restrictInputs") ) {
    restrictInputs_ = iConfig.getParameter<bool>("restrictInputs");
    maxInputs_      = iConfig.getParameter<unsigned int>("maxInputs");
  }
 

  string alias=iConfig.getUntrackedParameter<string>("alias",moduleLabel_);

  // make the "produces" statements
  makeProduces( alias, jetCollInstanceName_ );

  doFastJetNonUniform_ = false;
  if(iConfig.exists("doFastJetNonUniform")) doFastJetNonUniform_ = iConfig.getParameter<bool>   ("doFastJetNonUniform");
  if(doFastJetNonUniform_){
    puCenters_ = iConfig.getParameter<std::vector<double> >("puCenters");
    puWidth_ = iConfig.getParameter<double>("puWidth");
  }
  produces<std::vector<double> >("rhos");
  produces<std::vector<double> >("sigmas");
  produces<double>("rho");
  produces<double>("sigma");
  
}
VirtualJetProducer::~VirtualJetProducer ( ) [virtual]

Definition at line 226 of file VirtualJetProducer.cc.

{
} 

Member Function Documentation

void VirtualJetProducer::copyConstituents ( const std::vector< fastjet::PseudoJet > &  fjConstituents,
reco::Jet jet 
) [protected, virtual]

Definition at line 371 of file VirtualJetProducer.cc.

References reco::CompositePtrCandidate::addDaughter(), i, and inputs_.

{
  for (unsigned int i=0;i<fjConstituents.size();++i)
    jet->addDaughter(inputs_[fjConstituents[i].user_index()]);
}
vector< reco::CandidatePtr > VirtualJetProducer::getConstituents ( const std::vector< fastjet::PseudoJet > &  fjConstituents) [protected, virtual]

Definition at line 381 of file VirtualJetProducer.cc.

References i, getHLTprescales::index, inputs_, and query::result.

Referenced by FastjetJetProducer::produceTrackJets(), cms::SubEventGenJetProducer::runAlgorithm(), and writeJets().

{
  vector<reco::CandidatePtr> result;
  for (unsigned int i=0;i<fjConstituents.size();i++) {
    int index = fjConstituents[i].user_index();
    reco::CandidatePtr candidate = inputs_[index];
    result.push_back(candidate);
  }
  return result;
}
void VirtualJetProducer::inputTowers ( ) [protected, virtual]

Reimplemented in cms::CompoundJetProducer, cms::SubEventGenJetProducer, and SubjetFilterJetProducer.

Definition at line 318 of file VirtualJetProducer.cc.

References doPVCorrection_, fjInputs_, edm::Ptr< T >::get(), i, collect_tpl::input, inputEMin_, inputEtMin_, inputs_, isAnomalousTower(), edm::detail::isnan(), jetTypeE, makeCaloJet(), maxInputs_, CaloTower::p4(), restrictInputs_, python::multivaluedict::sort(), and vertex_.

Referenced by produce(), and FastjetJetProducer::produceTrackJets().

{
  std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
    inEnd = inputs_.end(), i = inBegin;
  for (; i != inEnd; ++i ) {
    reco::CandidatePtr input = *i;
    if (isnan(input->pt()))           continue;
    if (input->et()    <inputEtMin_)  continue;
    if (input->energy()<inputEMin_)   continue;
    if (isAnomalousTower(input))      continue;
    if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
      const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
      math::PtEtaPhiMLorentzVector ct(tower->p4(vertex_));
      fjInputs_.push_back(fastjet::PseudoJet(ct.px(),ct.py(),ct.pz(),ct.energy()));
    }
    else {
      fjInputs_.push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
                                            input->energy()));
    }
    fjInputs_.back().set_user_index(i - inBegin);
  }

  if ( restrictInputs_ && fjInputs_.size() > maxInputs_ ) {
    reco::helper::GreaterByPtPseudoJet   pTComparator;
    std::sort(fjInputs_.begin(), fjInputs_.end(), pTComparator);
    fjInputs_.resize(maxInputs_);
    edm::LogWarning("JetRecoTooManyEntries") << "Too many inputs in the event, limiting to first " << maxInputs_ << ". Output is suspect.";
  }
}
bool VirtualJetProducer::isAnomalousTower ( reco::CandidatePtr  input) [protected, virtual]
std::string VirtualJetProducer::jetType ( ) const [inline]

Definition at line 84 of file VirtualJetProducer.h.

References jetType_.

{ return jetType_; }
bool VirtualJetProducer::makeBasicJet ( const JetType::Type fTag) [inline, protected]

Definition at line 60 of file VirtualJetProducer.h.

References VirtualJetProducer::JetType::BasicJet.

                                                    {
    return fTag == JetType::BasicJet;
  }
bool VirtualJetProducer::makeCaloJet ( const JetType::Type fTag) [inline, protected]

Definition at line 48 of file VirtualJetProducer.h.

References VirtualJetProducer::JetType::CaloJet.

Referenced by inputTowers(), isAnomalousTower(), and produce().

                                                   {
    return fTag == JetType::CaloJet;
  }
bool VirtualJetProducer::makeGenJet ( const JetType::Type fTag) [inline, protected]

Definition at line 54 of file VirtualJetProducer.h.

References VirtualJetProducer::JetType::GenJet.

                                                  {
    return fTag == JetType::GenJet;
  }
bool VirtualJetProducer::makePFJet ( const JetType::Type fTag) [inline, protected]

Definition at line 51 of file VirtualJetProducer.h.

References VirtualJetProducer::JetType::PFJet.

                                                 {
    return fTag == JetType::PFJet;
  }
void VirtualJetProducer::makeProduces ( std::string  s,
std::string  tag = "" 
) [protected, virtual]

Definition at line 75 of file VirtualJetProducer.cc.

References GlobalPosition_Frontier_DevDB_cff::tag.

Referenced by SubjetFilterJetProducer::SubjetFilterJetProducer(), and VirtualJetProducer().

{
  if (makeCaloJet(jetTypeE)) {
    produces<reco::CaloJetCollection>(tag).setBranchAlias(alias);
  }
  else if (makePFJet(jetTypeE)) {
    produces<reco::PFJetCollection>(tag).setBranchAlias(alias);
  }
  else if (makeGenJet(jetTypeE)) {
    produces<reco::GenJetCollection>(tag).setBranchAlias(alias);
  }
  else if (makeTrackJet(jetTypeE)) {
    produces<reco::TrackJetCollection>(tag).setBranchAlias(alias);
  }
  else if (makeBasicJet(jetTypeE)) {
    produces<reco::BasicJetCollection>(tag).setBranchAlias(alias);
  }
}
bool VirtualJetProducer::makeTrackJet ( const JetType::Type fTag) [inline, protected]

Definition at line 57 of file VirtualJetProducer.h.

References VirtualJetProducer::JetType::TrackJet.

Referenced by FastjetJetProducer::produce().

                                                    {
    return fTag == JetType::TrackJet;
  }
void VirtualJetProducer::offsetCorrectJets ( std::vector< fastjet::PseudoJet > &  orphanInput) [protected]
void VirtualJetProducer::output ( edm::Event iEvent,
edm::EventSetup const &  iSetup 
) [protected, virtual]

Reimplemented in cms::CompoundJetProducer, and SubjetFilterJetProducer.

Definition at line 392 of file VirtualJetProducer.cc.

References VirtualJetProducer::JetType::BasicJet, VirtualJetProducer::JetType::CaloJet, Exception, VirtualJetProducer::JetType::GenJet, iEvent, jetTypeE, VirtualJetProducer::JetType::PFJet, and VirtualJetProducer::JetType::TrackJet.

Referenced by produce().

{
  // Write jets and constitutents. Will use fjJets_, inputs_
  // and fjClusterSeq_
  switch( jetTypeE ) {
  case JetType::CaloJet :
    writeJets<reco::CaloJet>( iEvent, iSetup);
    break;
  case JetType::PFJet :
    writeJets<reco::PFJet>( iEvent, iSetup);
    break;
  case JetType::GenJet :
    writeJets<reco::GenJet>( iEvent, iSetup);
    break;
  case JetType::TrackJet :
    writeJets<reco::TrackJet>( iEvent, iSetup);
    break;
  case JetType::BasicJet :
    writeJets<reco::BasicJet>( iEvent, iSetup);
    break;
  default:
    throw cms::Exception("InvalidInput") << "invalid jet type in VirtualJetProducer\n";
    break;
  };
  
}
void VirtualJetProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDProducer.

Reimplemented in cms::CATopJetProducer, CMSInsideOutJetProducer, FastjetJetProducer, cms::SubEventGenJetProducer, SubjetFilterJetProducer, and cms::SubJetProducer.

Definition at line 236 of file VirtualJetProducer.cc.

References doPUOffsetCorr_, doPVCorrection_, fjClusterSeq_, fjInputs_, fjJets_, edm::Event::getByLabel(), i, inputs_, inputTowers(), jetTypeE, LogDebug, makeCaloJet(), output(), runAlgorithm(), src_, srcPVs_, subtractor_, and vertex_.

{
  LogDebug("VirtualJetProducer") << "Entered produce\n";
  //determine signal vertex
  vertex_=reco::Jet::Point(0,0,0);
  if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
    LogDebug("VirtualJetProducer") << "Adding PV info\n";
    edm::Handle<reco::VertexCollection> pvCollection;
    iEvent.getByLabel(srcPVs_,pvCollection);
    if (pvCollection->size()>0) vertex_=pvCollection->begin()->position();
  }

  // For Pileup subtraction using offset correction:
  // set up geometry map
  if ( doPUOffsetCorr_ ) {
     subtractor_->setupGeometryMap(iEvent, iSetup);
  }

  // clear data
  LogDebug("VirtualJetProducer") << "Clear data\n";
  fjInputs_.clear();
  fjJets_.clear();
  inputs_.clear();  
  
  // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
  edm::Handle<reco::CandidateView> inputsHandle;
  iEvent.getByLabel(src_,inputsHandle);
  for (size_t i = 0; i < inputsHandle->size(); ++i) {
    inputs_.push_back(inputsHandle->ptrAt(i));
  }
  LogDebug("VirtualJetProducer") << "Got inputs\n";
  
  // Convert candidates to fastjet::PseudoJets.
  // Also correct to Primary Vertex. Will modify fjInputs_
  // and use inputs_
  fjInputs_.reserve(inputs_.size());
  inputTowers();
  LogDebug("VirtualJetProducer") << "Inputted towers\n";

  // For Pileup subtraction using offset correction:
  // Subtract pedestal. 
  if ( doPUOffsetCorr_ ) {
     subtractor_->reset(inputs_,fjInputs_,fjJets_);
     subtractor_->calculatePedestal(fjInputs_); 
     subtractor_->subtractPedestal(fjInputs_);    
     LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
  }

  // Run algorithm. Will modify fjJets_ and allocate fjClusterSeq_. 
  // This will use fjInputs_
  runAlgorithm( iEvent, iSetup );
  if ( doPUOffsetCorr_ ) {
     subtractor_->setAlgorithm(fjClusterSeq_);
  }

  LogDebug("VirtualJetProducer") << "Ran algorithm\n";

  // For Pileup subtraction using offset correction:
  // Now we find jets and need to recalculate their energy,
  // mark towers participated in jet,
  // remove occupied towers from the list and recalculate mean and sigma
  // put the initial towers collection to the jet,   
  // and subtract from initial towers in jet recalculated mean and sigma of towers 
  if ( doPUOffsetCorr_ ) {
    vector<fastjet::PseudoJet> orphanInput;
    subtractor_->calculateOrphanInput(orphanInput);
    subtractor_->calculatePedestal(orphanInput);
    subtractor_->offsetCorrectJets();
  }

  // Write the output jets.
  // This will (by default) call the member function template
  // "writeJets", but can be overridden. 
  // this will use inputs_
  output( iEvent, iSetup );
  LogDebug("VirtualJetProducer") << "Wrote jets\n";
  
  return;
}
virtual void VirtualJetProducer::runAlgorithm ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [protected, pure virtual]
template<typename T >
void VirtualJetProducer::writeJets ( edm::Event iEvent,
edm::EventSetup const &  iSetup 
) [protected]

Definition at line 420 of file VirtualJetProducer.cc.

References doAreaFastjet_, doFastJetNonUniform_, doPUOffsetCorr_, doRhoFastjet_, eta(), fjClusterSeq_, fjJets_, fjRangeDef_, getConstituents(), metsig::jet, analyzePatCleaning_cfg::jets, ExpressReco_HICollisions_FallBack::nEta, puCenters_, edm::Event::put(), puWidth_, dt_offlineAnalysis_common_cff::reco, rho, ExpressReco_HICollisions_FallBack::sigma, subtractor_, vertex_, and reco::writeSpecific().

{
  // produce output jet collection

  using namespace reco;

  std::auto_ptr<std::vector<T> > jets(new std::vector<T>() );
  jets->reserve(fjJets_.size());
      
  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
    // allocate this jet
    T jet;
    // get the fastjet jet
    const fastjet::PseudoJet& fjJet = fjJets_[ijet];
    // get the constituents from fastjet
    std::vector<fastjet::PseudoJet> fjConstituents =
      sorted_by_pt(fjClusterSeq_->constituents(fjJet));
    // convert them to CandidatePtr vector
    std::vector<CandidatePtr> constituents =
      getConstituents(fjConstituents);

    // calcuate the jet area
    double jetArea=0.0;
    if ( doAreaFastjet_ ) {
      fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
        dynamic_cast<fastjet::ClusterSequenceArea const *>(&*fjClusterSeq_);
      jetArea = clusterSequenceWithArea->area(fjJet);
    }
    
    // write the specifics to the jet (simultaneously sets 4-vector, vertex).
    // These are overridden functions that will call the appropriate
    // specific allocator. 
    writeSpecific(jet,
                  Particle::LorentzVector(fjJet.px(),
                                          fjJet.py(),
                                          fjJet.pz(),
                                          fjJet.E()),
                  vertex_, 
                  constituents, iSetup);
    
    jet.setJetArea (jetArea);

    if(doPUOffsetCorr_){
       jet.setPileup(subtractor_->getPileUpEnergy(ijet));
   }else{
       jet.setPileup (0.0);
    }

    // add to the list
    jets->push_back(jet);       
  }
  
  // put the jets in the collection
  iEvent.put(jets);
  
  if (doRhoFastjet_) {
    if(doFastJetNonUniform_){
      std::auto_ptr<std::vector<double> > rhos(new std::vector<double>);
      std::auto_ptr<std::vector<double> > sigmas(new std::vector<double>);
      int nEta = puCenters_.size();
      rhos->reserve(nEta);
      sigmas->reserve(nEta);
      fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
        dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );

      for(int ie = 0; ie < nEta; ++ie){
        double eta = puCenters_[ie];
        double rho = 0;
        double sigma = 0;
        double etamin=eta-puWidth_;
        double etamax=eta+puWidth_;
        fastjet::RangeDefinition range_rho(etamin,etamax);
        clusterSequenceWithArea->get_median_rho_and_sigma(range_rho, true, rho, sigma);
        rhos->push_back(rho);
        sigmas->push_back(sigma);
      }
      iEvent.put(rhos,"rhos");
      iEvent.put(sigmas,"sigmas");
    }else{
      std::auto_ptr<double> rho(new double(0.0));
      std::auto_ptr<double> sigma(new double(0.0));
      double mean_area = 0;
      
      fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
        dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );
      clusterSequenceWithArea->get_median_rho_and_sigma(*fjRangeDef_,false,*rho,*sigma,mean_area);
      iEvent.put(rho,"rho");
      iEvent.put(sigma,"sigma");
    }
    
  }
}

Member Data Documentation

Definition at line 159 of file VirtualJetProducer.h.

Referenced by VirtualJetProducer(), and writeJets().

Definition at line 162 of file VirtualJetProducer.h.

Referenced by produce(), VirtualJetProducer(), and writeJets().

Definition at line 149 of file VirtualJetProducer.h.

Referenced by inputTowers(), and produce().

Definition at line 178 of file VirtualJetProducer.h.

Referenced by FastjetJetProducer::runAlgorithm(), and VirtualJetProducer().

std::vector<fastjet::PseudoJet> VirtualJetProducer::fjInputs_ [protected]
std::vector<fastjet::PseudoJet> VirtualJetProducer::fjJets_ [protected]

Definition at line 177 of file VirtualJetProducer.h.

Referenced by VirtualJetProducer().

Definition at line 179 of file VirtualJetProducer.h.

Referenced by VirtualJetProducer(), and writeJets().

double VirtualJetProducer::inputEMin_ [protected]

Definition at line 147 of file VirtualJetProducer.h.

Referenced by inputTowers(), and cms::SubEventGenJetProducer::inputTowers().

double VirtualJetProducer::inputEtMin_ [protected]

Definition at line 146 of file VirtualJetProducer.h.

Referenced by inputTowers(), and cms::SubEventGenJetProducer::inputTowers().

std::string VirtualJetProducer::jetAlgorithm_ [protected]

Definition at line 144 of file VirtualJetProducer.h.

Referenced by VirtualJetProducer().

double VirtualJetProducer::jetPtMin_ [protected]
std::string VirtualJetProducer::jetType_ [protected]

Definition at line 143 of file VirtualJetProducer.h.

Referenced by jetType(), and VirtualJetProducer().

unsigned int VirtualJetProducer::maxBadEcalCells_ [protected]

Definition at line 166 of file VirtualJetProducer.h.

Referenced by isAnomalousTower().

unsigned int VirtualJetProducer::maxBadHcalCells_ [protected]

Definition at line 169 of file VirtualJetProducer.h.

Referenced by isAnomalousTower().

unsigned int VirtualJetProducer::maxInputs_ [protected]

Definition at line 153 of file VirtualJetProducer.h.

Referenced by inputTowers(), and VirtualJetProducer().

Definition at line 168 of file VirtualJetProducer.h.

Referenced by isAnomalousTower().

Definition at line 171 of file VirtualJetProducer.h.

Referenced by isAnomalousTower().

Definition at line 167 of file VirtualJetProducer.h.

Referenced by isAnomalousTower().

Definition at line 170 of file VirtualJetProducer.h.

Referenced by isAnomalousTower().

std::string VirtualJetProducer::moduleLabel_ [protected]
std::vector<double> VirtualJetProducer::puCenters_ [protected]

Definition at line 184 of file VirtualJetProducer.h.

Referenced by VirtualJetProducer(), and writeJets().

std::string VirtualJetProducer::puSubtractorName_ [protected]

Definition at line 163 of file VirtualJetProducer.h.

Referenced by VirtualJetProducer().

double VirtualJetProducer::puWidth_ [protected]

Definition at line 185 of file VirtualJetProducer.h.

Referenced by VirtualJetProducer(), and writeJets().

Definition at line 152 of file VirtualJetProducer.h.

Referenced by inputTowers(), and VirtualJetProducer().

double VirtualJetProducer::rParam_ [protected]

Definition at line 145 of file VirtualJetProducer.h.

Referenced by VirtualJetProducer().

Definition at line 142 of file VirtualJetProducer.h.

Referenced by produce(), and FastjetJetProducer::produceTrackJets().

boost::shared_ptr<PileUpSubtractor> VirtualJetProducer::subtractor_ [protected]

Definition at line 188 of file VirtualJetProducer.h.

Referenced by produce(), VirtualJetProducer(), and writeJets().