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

Added particle flow candidates to particle flow validation #29628

Merged
merged 16 commits into from May 5, 2020
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
17 changes: 8 additions & 9 deletions DQMOffline/Configuration/python/autoDQM.py
Expand Up @@ -116,7 +116,7 @@
'ctpps': ['DQMOfflineCTPPS',
'PostDQMOffline',
'DQMHarvestCTPPS'],

'btag': ['DQMOfflineBTag',
'PostDQMOffline',
'DQMHarvestBTag'],
Expand All @@ -141,8 +141,8 @@
'PostDQMOffline',
'HLTMonitoringClient'],

'HLTMonPA': ['HLTMonitoringPA',
'PostDQMOffline',
'HLTMonPA': ['HLTMonitoringPA',
'PostDQMOffline',
'HLTMonitoringClientPA'],

'express': ['@commonSiStripZeroBias+@muon+@hcal+@jetmet+@ecal',
Expand Down Expand Up @@ -177,12 +177,12 @@
'PostDQMOffline',
'DQMHarvestNanoAOD'],

'pfDQM': ['DQMOfflinePF',
'pfDQM': ['DQMOfflinePF+DQMOfflinePFExtended',
'PostDQMOffline',
'DQMHarvestPF'],

# 'standardDQM': ['@dcs+@DQMMessageLogger+@ecal+@hcal+@hcal2+@strip+@pixel+@castor+@ctpps+@muon+@tracking+@jetmet+@egamma+@L1TMon+@hlt+@btag+@beam+@physics+@HLTMon',
'standardDQM': ['DQMOffline',
# 'standardDQM': ['@dcs+@DQMMessageLogger+@ecal+@hcal+@hcal2+@strip+@pixel+@castor+@ctpps+@muon+@tracking+@jetmet+@egamma+@L1TMon+@hlt+@btag+@beam+@physics+@HLTMon',
'standardDQM': ['DQMOffline',
'PostDQMOffline',
'dqmHarvesting'],

Expand All @@ -194,8 +194,8 @@
'PostDQMOffline',
'dqmHarvestingExtraHLT'],

# 'standardDQMFakeHLT': ['@dcs+@DQMMessageLogger+@ecal+@hcal+@hcal2+@strip+@pixel+@castor+@ctpps+@muon+@tracking+@jetmet+@egamma+@L1TMon+@btag+@beam+@physics',
'standardDQMFakeHLT': ['DQMOfflineFakeHLT',
# 'standardDQMFakeHLT': ['@dcs+@DQMMessageLogger+@ecal+@hcal+@hcal2+@strip+@pixel+@castor+@ctpps+@muon+@tracking+@jetmet+@egamma+@L1TMon+@btag+@beam+@physics',
'standardDQMFakeHLT': ['DQMOfflineFakeHLT',
'PostDQMOffline',
'dqmHarvestingFakeHLT'],

Expand All @@ -217,4 +217,3 @@
for i in [0,2]:
autoDQM['phase2'][i] = '+'.join([autoDQM[m][i] for m in _phase2_allowed])
autoDQM['phase2'][1] = 'PostDQMOffline'

149 changes: 149 additions & 0 deletions Validation/RecoParticleFlow/plugins/PFCandidateAnalyzerDQM.cc
@@ -0,0 +1,149 @@
// Producer for particle flow candidates. Plots Eta, Phi, Charge, Pt (log freq, bin)
// for different types of particles described in python/defaults_cfi.py
// It actually uses packedCandidates so that we need only MINIAOD contents to run this DQMAnalyzer.
// note: for pt, log freq is done in this producer, but log freq is done by running
// compare.py
// author: Chosila Sutantawibul, April 23, 2020

#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"

#include "DataFormats/Common/interface/Handle.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/Event.h"

#include "DQMServices/Core/interface/DQMStore.h"
#include "DataFormats/PatCandidates/interface/PackedCandidate.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "DQMServices/Core/interface/DQMEDAnalyzer.h"

#include "TH1F.h"

#include <algorithm>
#include <cmath>
#include <fstream>
#include <iostream>
#include <memory>
#include <ostream>
#include <map>
#include <string>
#include <cstring>

class PFCandidateAnalyzerDQM : public DQMEDAnalyzer {
public:
explicit PFCandidateAnalyzerDQM(const edm::ParameterSet&);
void analyze(const edm::Event&, const edm::EventSetup&) override;

protected:
void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;

private:
//from config file
edm::InputTag PFCandTag;
edm::EDGetTokenT<edm::View<pat::PackedCandidate>> PFCandToken;
std::vector<double> etabins;
std::map<std::string, MonitorElement*> me;

std::map<uint32_t, std::string> pdgMap;
};

// constructor
PFCandidateAnalyzerDQM::PFCandidateAnalyzerDQM(const edm::ParameterSet& iConfig) {
PFCandTag = iConfig.getParameter<edm::InputTag>("PFCandType");
PFCandToken = consumes<edm::View<pat::PackedCandidate>>(PFCandTag);
etabins = iConfig.getParameter<std::vector<double>>("etabins");

//create map of pdgId
std::vector<uint32_t> pdgKeys = iConfig.getParameter<std::vector<uint32_t>>("pdgKeys");
std::vector<std::string> pdgStrs = iConfig.getParameter<std::vector<std::string>>("pdgStrs");
for (int i = 0, n = pdgKeys.size(); i < n; i++)
pdgMap[pdgKeys[i]] = pdgStrs[i];
}

void PFCandidateAnalyzerDQM::bookHistograms(DQMStore::IBooker& booker, edm::Run const&, edm::EventSetup const&) {
// all candidate
booker.setCurrentFolder("ParticleFlow/PackedCandidates/AllCandidate");

// for eta binning
int n = etabins.size() - 1;
float etabinArray[etabins.size()];
std::copy(etabins.begin(), etabins.end(), etabinArray);

//eta has variable bin sizes, use 4th def of TH1F constructor
TH1F* etaHist = new TH1F("AllCandidateEta", "AllCandidateEta", n, etabinArray);
me["AllCandidateEta"] = booker.book1D("AllCandidateEta", etaHist);

me["AllCandidateLog10Pt"] = booker.book1D("AllCandidateLog10Pt", "AllCandidateLog10Pt", 120, -2, 4);

//for phi binnings
double nPhiBins = 73;
double phiBinWidth = M_PI / (nPhiBins - 1) * 2.;
me["AllCandidatePhi"] = booker.book1D(
"AllCandidatePhi", "AllCandidatePhi", nPhiBins, -M_PI - 0.25 * phiBinWidth, +M_PI + 0.75 * phiBinWidth);

me["AllCandidateCharge"] = booker.book1D("AllCandidateCharge", "AllCandidateCharge", 3, -1.5, 1.5);
me["AllCandidatePtLow"] = booker.book1D("AllCandidatePtLow", "AllCandidatePtLow", 100, 0., 5.);
me["AllCandidatePtMid"] = booker.book1D("AllCandidatePtMid", "AllCandidatePtMid", 100, 0., 200.);
me["AllCandidatePtHigh"] = booker.book1D("AllCandidatePtHigh", "AllCandidatePtHigh", 100, 0., 1000.);

for (auto& pair : pdgMap) {
booker.setCurrentFolder("ParticleFlow/PackedCandidates/" + pair.second);

//TH1F only takes char*, so have to do conversions for histogram name
const char* etaHistName = (pair.second + "Eta").c_str();
TH1F* etaHist = new TH1F(etaHistName, etaHistName, n, etabinArray);
me[pair.second + "Eta"] = booker.book1D(pair.second + "Eta", etaHist);

me[pair.second + "Log10Pt"] = booker.book1D(pair.second + "Log10Pt", pair.second + "Log10Pt", 120, -2, 4);
me[pair.second + "Phi"] = booker.book1D(
pair.second + "Phi", pair.second + "Phi", nPhiBins, -M_PI - 0.25 * phiBinWidth, +M_PI + 0.75 * phiBinWidth);
me[pair.second + "Charge"] = booker.book1D(pair.second + "Charge", pair.second + "Charge", 3, -1.5, 1.5);
me[pair.second + "PtLow"] = booker.book1D(pair.second + "PtLow", pair.second + "PtLow", 100, 0., 5.);
me[pair.second + "PtMid"] = booker.book1D(pair.second + "PtMid", pair.second + "PtMid", 100, 0., 200.);
me[pair.second + "PtHigh"] = booker.book1D(pair.second + "PtHigh", pair.second + "PtHigh", 100, 0., 1000.);
}
}

void PFCandidateAnalyzerDQM::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
//retrieve
edm::Handle<edm::View<pat::PackedCandidate>> pfHandle;
iEvent.getByToken(PFCandToken, pfHandle);

if (!pfHandle.isValid()) {
edm::LogInfo("OutputInfo") << " failed to retrieve data required by ParticleFlow Task";
edm::LogInfo("OutputInfo") << " ParticleFlow Task cannot continue...!";
return;
} else {
//Analyze
// Loop Over Particle Flow Candidates

for (unsigned int i = 0; i < pfHandle->size(); i++) {
// Fill Histograms for Candidate Methods
// all candidates
me["AllCandidateLog10Pt"]->Fill(log10(pfHandle->at(i).pt()));
me["AllCandidateEta"]->Fill(pfHandle->at(i).eta());
me["AllCandidatePhi"]->Fill(pfHandle->at(i).phi());
me["AllCandidateCharge"]->Fill(pfHandle->at(i).charge());
me["AllCandidatePtLow"]->Fill(pfHandle->at(i).pt());
me["AllCandidatePtMid"]->Fill(pfHandle->at(i).pt());
me["AllCandidatePtHigh"]->Fill(pfHandle->at(i).pt());

int pdgId = abs(pfHandle->at(i).pdgId());
if (pdgMap.find(pdgId) != pdgMap.end()) {
me[pdgMap[pdgId] + "Log10Pt"]->Fill(log10(pfHandle->at(i).pt()));
me[pdgMap[pdgId] + "Eta"]->Fill(pfHandle->at(i).eta());
me[pdgMap[pdgId] + "Phi"]->Fill(pfHandle->at(i).phi());
me[pdgMap[pdgId] + "Charge"]->Fill(pfHandle->at(i).charge());
me[pdgMap[pdgId] + "PtLow"]->Fill(pfHandle->at(i).pt());
me[pdgMap[pdgId] + "PtMid"]->Fill(pfHandle->at(i).pt());
me[pdgMap[pdgId] + "PtHigh"]->Fill(pfHandle->at(i).pt());
}
}
}
}

#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(PFCandidateAnalyzerDQM);
13 changes: 12 additions & 1 deletion Validation/RecoParticleFlow/python/DQMForPF_MiniAOD_cff.py
Expand Up @@ -3,16 +3,27 @@
from Validation.RecoParticleFlow.particleFlowDQM_cff import pfJetAnalyzerDQM
from Validation.RecoParticleFlow.particleFlowDQM_cff import pfPuppiJetAnalyzerDQM
from Validation.RecoParticleFlow.particleFlowDQM_cff import pfJetDQMPostProcessor
from Validation.RecoParticleFlow.particleFlowDQM_cff import PFCandAnalyzerDQM
from Validation.RecoParticleFlow.offsetAnalyzerDQM_cff import offsetAnalyzerDQM
from Validation.RecoParticleFlow.offsetAnalyzerDQM_cff import offsetDQMPostProcessor
# Use also other POGs' analyzers for extended checks
from Validation.RecoMET.METRelValForDQM_cff import *
from Validation.RecoJets.JetValidation_cff import *

DQMOfflinePF = cms.Sequence(
pfJetAnalyzerDQM +
pfPuppiJetAnalyzerDQM +
offsetAnalyzerDQM
offsetAnalyzerDQM +
PFCandAnalyzerDQM
)

DQMHarvestPF = cms.Sequence(
pfJetDQMPostProcessor +
offsetDQMPostProcessor
)

# MET & Jets sequence
DQMOfflinePFExtended = cms.Sequence(
METValidationMiniAOD +
JetValidationMiniAOD
)
8 changes: 5 additions & 3 deletions Validation/RecoParticleFlow/python/defaults_cfi.py
Expand Up @@ -41,8 +41,8 @@ def genjet_distribution_name(ietabin):
Defined in PFCandidate.cc
ParticleTypes (use absolute value)
211 h, // charged hadron
11 e, // electron
13 mu, // muon
11 e, // electron
13 mu, // muon
22 gamma, // photon
130 h0, // neutral hadron
1 h_HF, // HF tower identified as a hadron
Expand All @@ -54,6 +54,9 @@ def genjet_distribution_name(ietabin):
#our particle types (note chu = unmatched charged hadrons)
candidateType = ["chm", "chu", "nh", "ne", "hfh", "hfe", "lep"]

# map reco::PFCandidate::pdgId() to particle type for PFCand
pdgIDDict = {211: "chargedHadron", 11:"electron", 13:"muon", 22:"photon", 130:"neutralHadron", 1:"HF_hadron", 2: "HF_EM_particle"}

offsetPlotBaseName = 'p_offset_eta'
offsetDir = 'ParticleFlow/Offset/'
offsetR = 0.4
Expand All @@ -62,4 +65,3 @@ def genjet_distribution_name(ietabin):

def offset_name( var, ivar, itype ) :
return "{0}_{1}{2}_{3}".format( offsetPlotBaseName, var, ivar, itype )

17 changes: 10 additions & 7 deletions Validation/RecoParticleFlow/python/particleFlowDQM_cff.py
@@ -1,4 +1,5 @@
import FWCore.ParameterSet.Config as cms
import Validation.RecoParticleFlow.defaults_cfi as default
from Validation.RecoParticleFlow.defaults_cfi import ptbins, etabins, response_distribution_name, genjet_distribution_name,jetResponseDir,genjetDir

#----- ----- ----- ----- ----- ----- ----- -----
Expand Down Expand Up @@ -94,12 +95,14 @@ def createGenJetPlots(ptbins, etabins):

)

#----- ----- ----- ----- ----- ----- ----- -----
#
# Sequence
# if you add things here, you also have to fix DQMForPF_MiniAOD_cff.py

pfDQM = cms.Sequence(
pfJetAnalyzerDQM *
pfPuppiJetAnalyzerDQM
# PFCandidates
PFCandAnalyzerDQM = cms.EDProducer("PFCandidateAnalyzerDQM",
PFCandType = cms.InputTag("packedPFCandidates"),
etabins = cms.vdouble( default.etaBinsOffset ),
pdgKeys = cms.vuint32( default.pdgIDDict.keys() ),
pdgStrs = cms.vstring( default.pdgIDDict.values() )
)


#----- ----- ----- ----- ----- ----- ----- -----