#include <RecoJets/FFTJetProducers/plugins/FFTJetCorrectionProducer.cc>
Public Member Functions | |
FFTJetCorrectionProducer (const edm::ParameterSet &) | |
~FFTJetCorrectionProducer () | |
Private Member Functions | |
template<typename Jet > | |
void | applyCorrections (edm::Event &iEvent, const edm::EventSetup &iSetup) |
virtual void | beginJob () |
virtual void | endJob () |
template<typename Jet > | |
void | makeProduces (const std::string &alias, const std::string &tag) |
template<typename Jet > | |
void | performPileupSubtraction (Jet &) |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
Private Attributes | |
unsigned long | eventCount |
const edm::InputTag | inputLabel |
const JetType | jetType |
const std::string | outputLabel |
const std::vector< std::string > | records |
std::vector< int > | sequenceMasks |
const bool | subtractPileup |
const bool | subtractPileupAs4Vec |
const bool | verbose |
const bool | writeUncertainties |
Description: producer for correcting jets created by FFTJetProducer
Implementation: [Notes on implementation]
Definition at line 93 of file FFTJetCorrectionProducer.cc.
FFTJetCorrectionProducer::FFTJetCorrectionProducer | ( | const edm::ParameterSet & | ps | ) | [explicit] |
Definition at line 361 of file FFTJetCorrectionProducer.cc.
References edm::ParameterSet::getUntrackedParameter(), jet_type_switch, makeProduces(), outputLabel, AlCaHLTBitMon_QueryRunRegistry::string, and writeUncertainties.
: inputLabel(ps.getParameter<edm::InputTag>("src")), outputLabel(ps.getParameter<std::string>("outputLabel")), jetType(parseJetType(ps.getParameter<std::string>("jetType"))), records(ps.getParameter<std::vector<std::string> >("records")), writeUncertainties(ps.getParameter<bool>("writeUncertainties")), subtractPileup(ps.getParameter<bool>("subtractPileup")), subtractPileupAs4Vec(ps.getParameter<bool>("subtractPileupAs4Vec")), verbose(ps.getUntrackedParameter<bool>("verbose", false)), eventCount(0UL) { const std::string alias(ps.getUntrackedParameter<std::string>( "alias", outputLabel)); jet_type_switch(makeProduces, alias, outputLabel); if (writeUncertainties) produces<std::vector<float> >(outputLabel).setBranchAlias(alias); }
FFTJetCorrectionProducer::~FFTJetCorrectionProducer | ( | ) |
Definition at line 381 of file FFTJetCorrectionProducer.cc.
{ }
void FFTJetCorrectionProducer::applyCorrections | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private] |
Definition at line 171 of file FFTJetCorrectionProducer.cc.
References coll, corr, gather_cfg::cout, Exception, f, reco::FFTJet< Real >::f_recoScale(), reco::FFTJet< Real >::f_status(), reco::FFTJet< Real >::f_vec(), edm::Event::getByLabel(), i, dtTPAnalyzer_cfg::inputLabel, instance, makeTFileFromDB::isMC, edm::EventBase::isRealData(), j, fwrapper::jets, testEve_cfg::level, PILEUP_CALCULATION_MASK, PILEUP_SUBTRACTION_MASK_ANY, edm::Event::put(), alignCSCRings::s, python::multivaluedict::sort(), and mathSSE::sqrt().
Referenced by produce().
{ using reco::FFTJet; // Various useful typedefs typedef reco::FFTAnyJet<Jet> MyJet; typedef std::vector<MyJet> MyCollection; typedef typename FFTJetCorrectorSequenceTypemap<MyJet>::loader Loader; typedef typename Loader::data_type CorrectorSequence; typedef typename CorrectorSequence::result_type CorrectorResult; // Load the jet corrector sequences const unsigned nRecords = records.size(); std::vector<edm::ESHandle<CorrectorSequence> > handles(nRecords); for (unsigned irec=0; irec<nRecords; ++irec) Loader::instance().load(iSetup, records[irec], handles[irec]); // Figure out which correction levels we are applying // and create masks which will indicate this sequenceMasks.clear(); sequenceMasks.reserve(nRecords); int totalMask = 0; for (unsigned irec=0; irec<nRecords; ++irec) { int levelMask = 0; const unsigned nLevels = handles[irec]->nLevels(); for (unsigned i=0; i<nLevels; ++i) { const unsigned lev = (*handles[irec])[i].level(); // Not tracking "level 0" corrections in the status word. // Level 0 is basically reserved for uncertainty calculations. if (lev) { const int mask = (1 << lev); if (totalMask & mask) throw cms::Exception("FFTJetBadConfig") << "Error in FFTJetCorrectionProducer::applyCorrections:" << " jet correction at level " << lev << " is applied more than once\n"; totalMask |= mask; levelMask |= mask; } } sequenceMasks.push_back(levelMask << 12); } totalMask = (totalMask << 12); // Is this data or MC? const bool isMC = !iEvent.isRealData(); // Load the jet collection edm::Handle<MyCollection> jets; iEvent.getByLabel(inputLabel, jets); // Create the output collection const unsigned nJets = jets->size(); std::auto_ptr<MyCollection> coll(new MyCollection()); coll->reserve(nJets); // Cycle over jets and apply the corrector sequences bool sorted = true; double previousPt = DBL_MAX; for (unsigned ijet=0; ijet<nJets; ++ijet) { const MyJet& j((*jets)[ijet]); // Check that this jet has not been corrected yet const int initialStatus = j.getFFTSpecific().f_status(); if (initialStatus & totalMask) throw cms::Exception("FFTJetBadConfig") << "Error in FFTJetCorrectionProducer::applyCorrections: " << "this jet collection is already corrected for some or all " << "of the specified levels\n"; MyJet corJ(j); if (verbose) { const reco::FFTJet<float>& fj = corJ.getFFTSpecific(); std::cout << "++++ Evt " << eventCount << " jet " << ijet << ": pt = " << corJ.pt() << ", eta = " << fj.f_vec().eta() << ", R = " << fj.f_recoScale() << ", s = 0x" << std::hex << fj.f_status() << std::dec << std::endl; } // Check if we need to subtract pileup first. // Pileup subtraction is not part of the corrector sequence // itself because 4-vector subtraction does not map well // into multiplication of 4-vectors by a scale factor. if (subtractPileup) { if (initialStatus & PILEUP_SUBTRACTION_MASK_ANY) throw cms::Exception("FFTJetBadConfig") << "Error in FFTJetCorrectionProducer::applyCorrections: " << "this jet collection is already pileup-subtracted\n"; if (!(initialStatus & PILEUP_CALCULATION_MASK)) throw cms::Exception("FFTJetBadConfig") << "Error in FFTJetCorrectionProducer::applyCorrections: " << "pileup was not calculated for this jet collection\n"; performPileupSubtraction(corJ); if (verbose) { const reco::FFTJet<float>& fj = corJ.getFFTSpecific(); std::cout << " Pileup subtract" << ": pt = " << corJ.pt() << ", eta = " << fj.f_vec().eta() << ", R = " << fj.f_recoScale() << ", s = 0x" << std::hex << fj.f_status() << std::dec << std::endl; } } // Apply all jet correction sequences double sigmaSquared = 0.0; for (unsigned irec=0; irec<nRecords; ++irec) { const CorrectorResult& corr = handles[irec]->correct(corJ, isMC); // Update the 4-vector FFTJet<float>& fftJet(const_cast<FFTJet<float>&>(corJ.getFFTSpecific())); corJ.setP4(corr.vec()); fftJet.setFourVec(corr.vec()); // Update the jet correction status fftJet.setStatus(fftJet.f_status() | sequenceMasks[irec]); // Update the (systematic) uncertainty const double s = corr.sigma(); sigmaSquared += s*s; } // There is no place for uncertainty in the jet structure. // However, there is the unused pileup field (FFTJet maintains // the pileup separately as a 4-vector). Use this unused field // to store the uncertainty. This hack is needed because // subsequent sequence sorting by Pt can change the jet ordering. if (writeUncertainties) corJ.setPileup(sqrt(sigmaSquared)); coll->push_back(corJ); // Check whether the sequence remains sorted by pt const double pt = corJ.pt(); if (pt > previousPt) sorted = false; previousPt = pt; if (verbose) { const reco::FFTJet<float>& fj = corJ.getFFTSpecific(); std::cout << " Fully corrected" << ": pt = " << corJ.pt() << ", eta = " << fj.f_vec().eta() << ", R = " << fj.f_recoScale() << ", s = 0x" << std::hex << fj.f_status() << std::dec << std::endl; } } if (!sorted) std::sort(coll->begin(), coll->end(), LocalSortByPt()); // Create the uncertainty sequence if (writeUncertainties) { std::auto_ptr<std::vector<float> > unc(new std::vector<float>()); unc->reserve(nJets); for (unsigned ijet=0; ijet<nJets; ++ijet) { MyJet& j((*coll)[ijet]); unc->push_back(j.pileup()); j.setPileup(0.f); } iEvent.put(unc, outputLabel); } iEvent.put(coll, outputLabel); ++eventCount; }
void FFTJetCorrectionProducer::beginJob | ( | void | ) | [private, virtual] |
void FFTJetCorrectionProducer::endJob | ( | void | ) | [private, virtual] |
void FFTJetCorrectionProducer::makeProduces | ( | const std::string & | alias, |
const std::string & | tag | ||
) | [private] |
Definition at line 144 of file FFTJetCorrectionProducer.cc.
References GlobalPosition_Frontier_DevDB_cff::tag.
Referenced by FFTJetCorrectionProducer().
{ produces<std::vector<reco::FFTAnyJet<T> > >(tag).setBranchAlias(alias); }
void FFTJetCorrectionProducer::performPileupSubtraction | ( | Jet & | jet | ) | [private] |
Definition at line 152 of file FFTJetCorrectionProducer.cc.
References fftjetcms::adjustForPileup(), PILEUP_SUBTRACTION_MASK_4VEC, PILEUP_SUBTRACTION_MASK_PT, and ntuplemaker::status.
{ using reco::FFTJet; FFTJet<float>& fftJet(const_cast<FFTJet<float>&>(jet.getFFTSpecific())); const math::XYZTLorentzVector& new4vec = adjustForPileup( fftJet.f_vec(), fftJet.f_pileup(), subtractPileupAs4Vec); fftJet.setFourVec(new4vec); int status = fftJet.f_status(); if (subtractPileupAs4Vec) status |= PILEUP_SUBTRACTION_MASK_4VEC; else status |= PILEUP_SUBTRACTION_MASK_PT; fftJet.setStatus(status); jet.setP4(new4vec); }
void FFTJetCorrectionProducer::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 387 of file FFTJetCorrectionProducer.cc.
References applyCorrections(), and jet_type_switch.
{ jet_type_switch(applyCorrections, iEvent, iSetup); }
unsigned long FFTJetCorrectionProducer::eventCount [private] |
Definition at line 139 of file FFTJetCorrectionProducer.cc.
const edm::InputTag FFTJetCorrectionProducer::inputLabel [private] |
Definition at line 114 of file FFTJetCorrectionProducer.cc.
const JetType FFTJetCorrectionProducer::jetType [private] |
Definition at line 120 of file FFTJetCorrectionProducer.cc.
const std::string FFTJetCorrectionProducer::outputLabel [private] |
Definition at line 117 of file FFTJetCorrectionProducer.cc.
Referenced by FFTJetCorrectionProducer().
const std::vector<std::string> FFTJetCorrectionProducer::records [private] |
Definition at line 123 of file FFTJetCorrectionProducer.cc.
std::vector<int> FFTJetCorrectionProducer::sequenceMasks [private] |
Definition at line 136 of file FFTJetCorrectionProducer.cc.
const bool FFTJetCorrectionProducer::subtractPileup [private] |
Definition at line 129 of file FFTJetCorrectionProducer.cc.
const bool FFTJetCorrectionProducer::subtractPileupAs4Vec [private] |
Definition at line 130 of file FFTJetCorrectionProducer.cc.
const bool FFTJetCorrectionProducer::verbose [private] |
Definition at line 133 of file FFTJetCorrectionProducer.cc.
const bool FFTJetCorrectionProducer::writeUncertainties [private] |
Definition at line 126 of file FFTJetCorrectionProducer.cc.
Referenced by FFTJetCorrectionProducer().