Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix calculation of smeared MET on MiniAOD #36165

Closed
18 changes: 18 additions & 0 deletions DataFormats/PatCandidates/interface/Jet.h
Original file line number Diff line number Diff line change
Expand Up @@ -177,6 +177,20 @@ namespace pat {
void scaleEnergy(double fScale) override { scaleEnergy(fScale, "Unscaled"); }
void scaleEnergy(double fScale, const std::string& level);

// Methods to work with a JER correction factor
//
// is the JER correction factor valid, i.e. has it been written using the saveJerFactor method
bool isJerFactorValid() const;

// load the currently saved JER correction factor and check whether the factor was propely set before loading
float jerFactor() const;

// save a JER correction factor
void setJerFactor(float jerFactor);

// reset a saved JER correction factor
void resetJerFactor();
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I'm not quite sure where this is needed later?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

At the moment this is not needed. I just added it for convenience. In the new commit, it is still included, however I can also remove it if desired.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I would just remove it to avoid confusion.

Copy link
Contributor Author

@michaelwassmer michaelwassmer Dec 8, 2021

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Hi again,
I now implemented it in a way to better keep track whether the jerFactor is already set or not. With this the user should be aware that he overwrites an already existing factor and can then decide whether this is correct or not. If so, the resetJerFactor method can be used to be able to overwrite. Hope this is better and the last commit.


private:
/// index of the set of jec factors with given label; returns -1 if no set
/// of jec factors exists with the given label
Expand Down Expand Up @@ -650,6 +664,10 @@ namespace pat {
// ---- id functions ----
reco::JetID jetID_;

// ---- JER correction factor ----
float jerFactor_;
bool jerFactorValid_;

private:
// ---- helper functions ----

Expand Down
61 changes: 57 additions & 4 deletions DataFormats/PatCandidates/src/Jet.cc
Original file line number Diff line number Diff line change
Expand Up @@ -10,23 +10,43 @@ using namespace pat;

/// default constructor
Jet::Jet()
: PATObject<reco::Jet>(reco::Jet()), embeddedCaloTowers_(false), embeddedPFCandidates_(false), jetCharge_(0.) {}
: PATObject<reco::Jet>(reco::Jet()),
embeddedCaloTowers_(false),
embeddedPFCandidates_(false),
jetCharge_(0.),
jerFactor_(0.0),
jerFactorValid_(false) {}

/// constructor from a reco::Jet
Jet::Jet(const reco::Jet& aJet)
: PATObject<reco::Jet>(aJet), embeddedCaloTowers_(false), embeddedPFCandidates_(false), jetCharge_(0.0) {
: PATObject<reco::Jet>(aJet),
embeddedCaloTowers_(false),
embeddedPFCandidates_(false),
jetCharge_(0.0),
jerFactor_(0.0),
jerFactorValid_(false) {
tryImportSpecific(aJet);
}

/// constructor from ref to reco::Jet
Jet::Jet(const edm::Ptr<reco::Jet>& aJetRef)
: PATObject<reco::Jet>(aJetRef), embeddedCaloTowers_(false), embeddedPFCandidates_(false), jetCharge_(0.0) {
: PATObject<reco::Jet>(aJetRef),
embeddedCaloTowers_(false),
embeddedPFCandidates_(false),
jetCharge_(0.0),
jerFactor_(0.0),
jerFactorValid_(false) {
tryImportSpecific(*aJetRef);
}

/// constructor from ref to reco::Jet
Jet::Jet(const edm::RefToBase<reco::Jet>& aJetRef)
: PATObject<reco::Jet>(aJetRef), embeddedCaloTowers_(false), embeddedPFCandidates_(false), jetCharge_(0.0) {
: PATObject<reco::Jet>(aJetRef),
embeddedCaloTowers_(false),
embeddedPFCandidates_(false),
jetCharge_(0.0),
jerFactor_(0.0),
jerFactorValid_(false) {
tryImportSpecific(*aJetRef);
}

Expand Down Expand Up @@ -596,3 +616,36 @@ void Jet::addSubjets(pat::JetPtrCollection const& pieces, std::string const& lab
subjetCollections_.push_back(pieces);
subjetLabels_.push_back(label);
}

// is the JER correction factor valid, i.e. has it been written using the saveJerFactor method
bool Jet::isJerFactorValid() const { return jerFactorValid_; }

// load the currently saved JER correction factor and check whether the factor was propely set before loading
float Jet::jerFactor() const {
if (isJerFactorValid()) {
return jerFactor_;
} else {
throw cms::Exception(
"Loading JER factor of a pat::Jet, but it was not set properly before. You need to save a JER factor first, "
"see methods in pat::Jet.");
return 0.0;
}
}

// save a JER correction factor
void Jet::setJerFactor(float jerFactor) {
if (isJerFactorValid()) {
throw cms::Exception(
"Setting a JER factor although the JER factor was already set before."
"Please reset the JER factor before setting a new one, see methods in pat::Jet");
} else {
jerFactor_ = jerFactor;
jerFactorValid_ = true;
}
}

// reset a saved JER correction factor
void Jet::resetJerFactor() {
jerFactor_ = 0.0;
jerFactorValid_ = false;
}
4 changes: 2 additions & 2 deletions DataFormats/PatCandidates/src/MET.cc
Original file line number Diff line number Diff line change
Expand Up @@ -342,9 +342,9 @@ void MET::setUncShift(double px, double py, double sumEt, METUncertainty shift,
//changing reference to only get the uncertainty shift and not the smeared one
// which is performed independently
shift = (MET::METUncertainty)(METUncertainty::METUncertaintySize + shift + 1);
const PackedMETUncertainty &ref = uncertainties_[METUncertainty::NoShift];
const PackedMETUncertainty &ref = corrections_[METCorrectionType::Smear];
uncertainties_[shift].set(
px + ref.dpx() - this->px(), py + ref.dpy() - this->py(), sumEt + ref.dsumEt() - this->sumEt());
px - ref.dpx() - this->px(), py - ref.dpy() - this->py(), sumEt - ref.dsumEt() - this->sumEt());
} else
uncertainties_[shift].set(px - this->px(), py - this->py(), sumEt - this->sumEt());
}
Expand Down
3 changes: 2 additions & 1 deletion DataFormats/PatCandidates/src/classes_def_objects.xml
Original file line number Diff line number Diff line change
Expand Up @@ -248,7 +248,8 @@
newObj->setPflowIsolationVariables(pfIsoVar);]]>
</ioread>

<class name="pat::Jet" ClassVersion="17">
<class name="pat::Jet" ClassVersion="18">
<version ClassVersion="18" checksum="3195476741"/>
<version ClassVersion="17" checksum="739868501"/>
<field name="daughtersTemp_" transient="true"/>
<version ClassVersion="16" checksum="4069285947"/>
Expand Down
14 changes: 14 additions & 0 deletions JetMETCorrections/Type1MET/interface/PFJetMETcorrInputProducerT.h
Original file line number Diff line number Diff line change
Expand Up @@ -61,6 +61,17 @@ namespace PFJetMETcorrInputProducer_namespace {
reco::Candidate::LorentzVector operator()(const T& jet) const { return jet.p4(); }
};

// functor to retrieve JER factors of jets
// for pat::Jet a previously saved JER factor is retrieved
// general template just returns 1
// specialized template for pat::Jet returns desired JER correction value if available
template <typename T>
class RetrieveJerT {
public:
RetrieveJerT() {}
float operator()(const T& jet) const { return 1.0; }
};

} // namespace PFJetMETcorrInputProducer_namespace

template <typename T, typename Textractor>
Expand Down Expand Up @@ -196,6 +207,9 @@ class PFJetMETcorrInputProducerT : public edm::stream::EDProducer<> {
corrJetP4 = jetCorrExtractor_(jet, jetCorrLabel_.label(), jetCorrEtaMax_, &rawJetP4);
else
corrJetP4 = jetCorrExtractor_(jet, jetCorr.product(), jetCorrEtaMax_, &rawJetP4);
// retrieve JER factors in case of pat::Jets (done via specialized template defined for pat::Jets) and apply it
const static PFJetMETcorrInputProducer_namespace::RetrieveJerT<T> retrieveJER{};
corrJetP4 *= retrieveJER(jet);

if (corrJetP4.pt() > type1JetPtThreshold_) {
reco::Candidate::LorentzVector rawJetP4offsetCorr = rawJetP4;
Expand Down
18 changes: 17 additions & 1 deletion PhysicsTools/PatUtils/interface/SmearedJetProducerT.h
Original file line number Diff line number Diff line change
Expand Up @@ -26,6 +26,7 @@
#include "DataFormats/JetReco/interface/GenJet.h"
#include "DataFormats/JetReco/interface/GenJetCollection.h"
#include "DataFormats/Math/interface/deltaR.h"
#include "DataFormats/PatCandidates/interface/Jet.h"

#include "JetMETCorrections/Modules/interface/JetResolution.h"

Expand Down Expand Up @@ -93,6 +94,21 @@ namespace pat {
};
}; // namespace pat

namespace SmearedJetProducer_namespace {
// template function to apply JER
template <typename T>
void SmearJet(T& jet, float smearfactor) {
jet.scaleEnergy(smearfactor);
}
// template specialization for pat::Jets to store the JER factor
template <>
inline void SmearJet<pat::Jet>(pat::Jet& jet, float smearfactor) {
jet.scaleEnergy(smearfactor);
jet.resetJerFactor();
jet.setJerFactor(smearfactor);
}
} // namespace SmearedJetProducer_namespace

template <typename T>
class SmearedJetProducerT : public edm::stream::EDProducer<> {
using JetCollection = std::vector<T>;
Expand Down Expand Up @@ -288,7 +304,7 @@ class SmearedJetProducerT : public edm::stream::EDProducer<> {
}

T smearedJet = jet;
smearedJet.scaleEnergy(smearFactor);
SmearedJetProducer_namespace::SmearJet(smearedJet, smearFactor);

if (m_debug) {
std::cout << "smeared jet (" << smearFactor << "): pt: " << smearedJet.pt() << " eta: " << smearedJet.eta()
Expand Down
28 changes: 26 additions & 2 deletions PhysicsTools/PatUtils/plugins/PATPFJetMETcorrInputProducer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -16,17 +16,41 @@ namespace PFJetMETcorrInputProducer_namespace {
bool isPatJet(const pat::Jet& jet) const { return true; }
};

// template specialization for pat::Jet
// remove JER correction if JER factor was saved previously
template <>
class RawJetExtractorT<pat::Jet> {
public:
RawJetExtractorT() {}
reco::Candidate::LorentzVector operator()(const pat::Jet& jet) const {
reco::Candidate::LorentzVector uncorrected_jet;
if (jet.jecSetsAvailable())
return jet.correctedP4("Uncorrected");
uncorrected_jet = jet.correctedP4("Uncorrected");
else
return jet.p4();
uncorrected_jet = jet.p4();
// remove JER correction factor from pat::Jets
if (jet.isJerFactorValid()) {
uncorrected_jet *= (1.0 / jet.jerFactor());
}
return uncorrected_jet;
}
};

// template specialization for pat::Jets
// retrieve JER factor if it was saved previously
// otherwise just return 1
template <>
class RetrieveJerT<pat::Jet> {
public:
RetrieveJerT() {}
float operator()(const pat::Jet& jet) const {
if (jet.isJerFactorValid()) {
return jet.jerFactor();
} else
return 1.0;
}
};

} // namespace PFJetMETcorrInputProducer_namespace

typedef PFJetMETcorrInputProducerT<pat::Jet, PATJetCorrExtractor> PATPFJetMETcorrInputProducer;
Expand Down