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

Update of MET significance tuning for Phys14 samples #8396

Merged
merged 2 commits into from Mar 19, 2015
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
6 changes: 2 additions & 4 deletions RecoMET/METAlgorithms/interface/METSignificance.h
Expand Up @@ -36,15 +36,13 @@ namespace metsig {
~METSignificance();

reco::METCovMatrix getCovariance(const edm::View<reco::Jet>& jets,
const std::vector<reco::Candidate::LorentzVector>& leptons,
const std::vector< edm::Handle<reco::CandidateView> >& leptons,
const edm::View<reco::Candidate>& pfCandidates);
double getSignificance(const reco::METCovMatrix& cov, const reco::MET& met ) const;

private:
std::vector<reco::Jet> cleanJets(const edm::View<reco::Jet>& jets,
const std::vector<reco::Candidate::LorentzVector>& leptons);
bool cleanJet(const reco::Jet& jet,
const std::vector<reco::Candidate::LorentzVector>& leptons);
const std::vector< edm::Handle<reco::CandidateView> >& leptons );

double jetThreshold_;
double dR2match_;
Expand Down
86 changes: 56 additions & 30 deletions RecoMET/METAlgorithms/src/METSignificance.cc
Expand Up @@ -50,31 +50,62 @@ metsig::METSignificance::~METSignificance() {

reco::METCovMatrix
metsig::METSignificance::getCovariance(const edm::View<reco::Jet>& jets,
const std::vector<reco::Candidate::LorentzVector>& leptons,
const std::vector< edm::Handle<reco::CandidateView> >& leptons,
const edm::View<reco::Candidate>& pfCandidates) {

// metsig covariance
double cov_xx = 0;
double cov_xy = 0;
double cov_yy = 0;

// for lepton and jet subtraction
std::vector<reco::CandidatePtr> footprint;

// subtract leptons out of sumPt
for ( std::vector< edm::Handle<reco::CandidateView> >::const_iterator lep_i = leptons.begin();
lep_i != leptons.end(); ++lep_i ) {
for( reco::CandidateView::const_iterator lep = (*lep_i)->begin(); lep != (*lep_i)->end(); lep++ ){
if( lep->pt() > 10 ){
for( unsigned int n=0; n < lep->numberOfSourceCandidatePtrs(); n++ ){
if( lep->sourceCandidatePtr(n).isNonnull() and lep->sourceCandidatePtr(n).isAvailable() ){
footprint.push_back(lep->sourceCandidatePtr(n));
}
}
}
}
}
// subtract jets out of sumPt
for(edm::View<reco::Jet>::const_iterator jet = jets.begin(); jet != jets.end(); ++jet) {

// disambiguate jets and leptons
if(!cleanJet(*jet, leptons) ) continue;

for( unsigned int n=0; n < jet->numberOfSourceCandidatePtrs(); n++){
if( jet->sourceCandidatePtr(n).isNonnull() and jet->sourceCandidatePtr(n).isAvailable() ){
footprint.push_back(jet->sourceCandidatePtr(n));
}
}

}

// calculate sumPt
double sumPt = 0;
for( edm::View<reco::Candidate>::const_iterator cand = pfCandidates.begin();
cand != pfCandidates.end(); ++cand){
sumPt += cand->pt();
}

//MM & NM :subtraction of lepton will need more study to see which option is optimal
//does not signficantly change the results
// check if candidate exists in a lepton or jet
bool cleancand = true;
for(unsigned int i=0; i < footprint.size(); i++){
if( footprint[i]->p4() == cand->p4() ){
cleancand = false;
}
}
// if not, add to sumPt
if( cleancand ){
sumPt += cand->pt();
}

// subtract leptons out of sumPt
// for ( std::vector<reco::Candidate::LorentzVector>::const_iterator lepton = leptons.begin();
// lepton != leptons.end(); ++lepton ) {
// //sumPt -= lepton->Pt();
// }
// //protection against unphysical events
// if(sumPt<0) sumPt=0;
}

// add jets to metsig covariance matrix and subtract them from sumPt
for(edm::View<reco::Jet>::const_iterator jet = jets.begin(); jet != jets.end(); ++jet) {
Expand Down Expand Up @@ -111,24 +142,19 @@ metsig::METSignificance::getCovariance(const edm::View<reco::Jet>& jets,
cov_xy += (dpt*dpt-dph*dph)*c*s;
cov_yy += dph*dph*c*c + dpt*dpt*s*s;

// subtract the pf constituents in each jet out of the sumPt
for(unsigned int i=0; i < jet->numberOfDaughters(); i++){
sumPt -= jet->daughter(i)->pt();
}

} else {

// subtract the pf constituents in each jet out of the sumPt
for(unsigned int i=0; i < jet->numberOfDaughters(); i++){
sumPt -= jet->daughter(i)->pt();
}
// add the (corrected) jet to the sumPt
sumPt += jpt;

}

}


//protection against unphysical events
if(sumPt<0) sumPt=0;

// add pseudo-jet to metsig covariance matrix
cov_xx += pjetParams_[0]*pjetParams_[0] + pjetParams_[1]*pjetParams_[1]*sumPt;
cov_yy += pjetParams_[0]*pjetParams_[0] + pjetParams_[1]*pjetParams_[1]*sumPt;
Expand Down Expand Up @@ -161,15 +187,15 @@ metsig::METSignificance::getSignificance(const reco::METCovMatrix& cov, const re

bool
metsig::METSignificance::cleanJet(const reco::Jet& jet,
const std::vector<reco::Candidate::LorentzVector>& leptons) {
if ( jet.pt() < jetThreshold_ ) return false;

for ( std::vector<reco::Candidate::LorentzVector>::const_iterator lepton = leptons.begin();
lepton != leptons.end(); ++lepton ) {
if ( reco::deltaR2( *lepton, jet) < dR2match_ ) {
return false;
}
const std::vector< edm::Handle<reco::CandidateView> >& leptons ){

for ( std::vector< edm::Handle<reco::CandidateView> >::const_iterator lep_i = leptons.begin();
lep_i != leptons.end(); ++lep_i ) {
for( reco::CandidateView::const_iterator lep = (*lep_i)->begin(); lep != (*lep_i)->end(); lep++ ){
if ( reco::deltaR2(*lep, jet) < dR2match_ ) {
return false;
}
}
}
return true;
}
4 changes: 2 additions & 2 deletions RecoMET/METProducers/python/METSignificanceParams_cfi.py
Expand Up @@ -14,7 +14,7 @@
jeta = cms.vdouble(0.5, 1.1, 1.7, 2.3),

# tuning parameters
jpar = cms.vdouble(1.15061, 1.07776, 1.04204, 1.12509, 1.56414),
pjpar = cms.vdouble(0.0, 0.548758)
jpar = cms.vdouble(1.41,1.29,1.41,1.40,2.53),
pjpar = cms.vdouble(0.0,0.674)

)
12 changes: 5 additions & 7 deletions RecoMET/METProducers/src/METSignificanceProducer.cc
Expand Up @@ -53,16 +53,14 @@ namespace cms
//
// leptons
//
std::vector<reco::Candidate::LorentzVector> leptons;
std::vector< edm::Handle<reco::CandidateView> > leptons;
for ( std::vector<edm::EDGetTokenT<edm::View<reco::Candidate> > >::const_iterator srcLeptons_i = lepTokens_.begin();
srcLeptons_i != lepTokens_.end(); ++srcLeptons_i ) {

edm::Handle<reco::CandidateView> leptons_i;
edm::Handle<reco::CandidateView> leptons_i;
event.getByToken(*srcLeptons_i, leptons_i);
for ( reco::CandidateView::const_iterator lepton = leptons_i->begin();
lepton != leptons_i->end(); ++lepton ) {
leptons.push_back(lepton->p4());
}
leptons.push_back( leptons_i );

}

//
Expand All @@ -87,7 +85,7 @@ namespace cms

event.put( covPtr, "METCovariance" );
event.put( significance, "METSignificance" );

}

//____________________________________________________________________________||
Expand Down
7 changes: 5 additions & 2 deletions RecoMET/METProducers/src/PFMETProducer.cc
Expand Up @@ -70,15 +70,18 @@ namespace cms
reco::METCovMatrix PFMETProducer::getMETCovMatrix(const edm::Event& event, const edm::View<reco::Candidate>& candInput) const {

// leptons
std::vector<reco::Candidate::LorentzVector> leptons;
std::vector< edm::Handle<reco::CandidateView> > leptons;
for ( std::vector<edm::EDGetTokenT<edm::View<reco::Candidate> > >::const_iterator srcLeptons_i = lepTokens_.begin();
srcLeptons_i != lepTokens_.end(); ++srcLeptons_i ) {
edm::Handle<reco::CandidateView> leptons_i;
event.getByToken(*srcLeptons_i, leptons_i);
leptons.push_back( leptons_i );
/*
for ( reco::CandidateView::const_iterator lepton = leptons_i->begin();
lepton != leptons_i->end(); ++lepton ) {
leptons.push_back(lepton->p4());
leptons.push_back(*lepton);
}
*/
}

// jets
Expand Down
5 changes: 3 additions & 2 deletions RecoMET/METProducers/test/recoMET_METSignificance_cfg.py
Expand Up @@ -16,7 +16,8 @@
"PoolSource",
fileNames = cms.untracked.vstring(
### MINIAODSIM
'/store/relval/CMSSW_7_3_0_pre1/RelValTTbar_13/MINIAODSIM/PU25ns_PRE_LS172_V15-v1/00000//A652D13A-335F-E411-90BA-02163E008D01.root'
'/store/relval/CMSSW_7_3_0/RelValZMM_13/MINIAODSIM/MCRUN2_73_V7-v1/00000/127CA68E-8981-E411-A524-002590593872.root',
'/store/relval/CMSSW_7_3_0/RelValZMM_13/MINIAODSIM/MCRUN2_73_V7-v1/00000/56FE228D-8981-E411-9AD8-0025905A6126.root'
)
)

Expand All @@ -34,7 +35,7 @@
##____________________________________________________________________________||
process.options = cms.untracked.PSet(wantSummary = cms.untracked.bool(True))
process.MessageLogger.cerr.FwkReport.reportEvery = 50
process.maxEvents = cms.untracked.PSet(input = cms.untracked.int32(-1))
process.maxEvents = cms.untracked.PSet(input = cms.untracked.int32(100))

##____________________________________________________________________________||
process.p = cms.Path(
Expand Down