-
Notifications
You must be signed in to change notification settings - Fork 4.4k
Commit
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Merge pull request #34384 from SissonJ/CMSSW_EfficiencyExtras_Updated
DQM Pixel Phase 1: Several new Efficiency trend plots
- Loading branch information
Showing
5 changed files
with
273 additions
and
5 deletions.
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
246 changes: 246 additions & 0 deletions
246
DQM/SiPixelPhase1Track/plugins/SiPixelPhase1EfficiencyExtras.cc
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,246 @@ | ||
// -*- C++ -*- | ||
// | ||
// Package: SiPixelPhase1EfficiencyExtras | ||
// Class: SiPixelPhase1EfficiencyExtras | ||
// | ||
/**\class | ||
Description: Create the Phase 1 extra efficiency trend plots | ||
Implementation: | ||
<Notes on implementation> | ||
*/ | ||
// | ||
// Original Author: Jack Sisson, Julie Hogan | ||
// Created: 7 July, 2021 | ||
// | ||
// | ||
|
||
// Framework | ||
#include "FWCore/ServiceRegistry/interface/Service.h" | ||
#include "FWCore/MessageLogger/interface/MessageLogger.h" | ||
#include "FWCore/ParameterSet/interface/ParameterSet.h" | ||
#include "FWCore/Framework/interface/Frameworkfwd.h" | ||
#include "FWCore/Framework/interface/MakerMacros.h" | ||
#include "FWCore/Framework/interface/Event.h" | ||
#include "FWCore/Framework/interface/EventSetup.h" | ||
#include "DataFormats/Common/interface/Handle.h" | ||
|
||
#include "FWCore/ParameterSet/interface/ParameterSet.h" | ||
|
||
// DQM Framework | ||
#include "DQMServices/Core/interface/DQMStore.h" | ||
#include "DQMServices/Core/interface/DQMEDHarvester.h" | ||
|
||
using namespace std; | ||
using namespace edm; | ||
|
||
class SiPixelPhase1EfficiencyExtras : public DQMEDHarvester { | ||
public: | ||
explicit SiPixelPhase1EfficiencyExtras(const edm::ParameterSet& conf); | ||
~SiPixelPhase1EfficiencyExtras() override; | ||
|
||
protected: | ||
void beginRun(edm::Run const& run, edm::EventSetup const& eSetup) override; | ||
|
||
void dqmEndJob(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter) override; | ||
|
||
std::string effFolderName_; | ||
std::string vtxFolderName_; | ||
std::string instLumiFolderName_; | ||
|
||
private: | ||
edm::ParameterSet conf_; | ||
}; | ||
|
||
SiPixelPhase1EfficiencyExtras::SiPixelPhase1EfficiencyExtras(const edm::ParameterSet& iConfig) : conf_(iConfig) { | ||
LogInfo("PixelDQM") << "SiPixelPhase1EfficiencyExtras::SiPixelPhase1EfficiencyExtras: Hello!" << endl; | ||
effFolderName_ = conf_.getParameter<std::string>("EffFolderName"); | ||
vtxFolderName_ = conf_.getParameter<std::string>("VtxFolderName"); | ||
instLumiFolderName_ = conf_.getParameter<std::string>("InstLumiFolderName"); | ||
} | ||
|
||
SiPixelPhase1EfficiencyExtras::~SiPixelPhase1EfficiencyExtras() { | ||
// do anything here that needs to be done at desctruction time | ||
// (e.g. close files, deallocate resources etc.) | ||
LogInfo("PixelDQM") << "SiPixelPhase1EfficiencyExtras::~SiPixelPhase1EfficiencyExtras: Destructor" << endl; | ||
} | ||
|
||
void SiPixelPhase1EfficiencyExtras::beginRun(edm::Run const& run, edm::EventSetup const& eSetup) {} | ||
|
||
//------------------------------------------------------------------ | ||
// Method called for every event | ||
//------------------------------------------------------------------ | ||
void SiPixelPhase1EfficiencyExtras::dqmEndJob(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter) { | ||
iBooker.setCurrentFolder(effFolderName_); | ||
|
||
//Get the existing histos | ||
MonitorElement* vtx_v_lumi = iGetter.get(vtxFolderName_ + "/NumberOfGoodPVtxVsLS_GenTk"); | ||
|
||
MonitorElement* scalLumi_v_lumi = iGetter.get(instLumiFolderName_ + "/lumiVsLS"); | ||
|
||
MonitorElement* eff_v_lumi_forward = | ||
iGetter.get(effFolderName_ + "/hitefficiency_per_Lumisection_per_PXDisk_PXForward"); | ||
|
||
MonitorElement* eff_v_lumi_barrel = | ||
iGetter.get(effFolderName_ + "/hitefficiency_per_Lumisection_per_PXLayer_PXBarrel"); | ||
|
||
//set up some booleans that will tell us which graphs to create | ||
bool createNvtx = true; | ||
bool createInstLumi = true; | ||
|
||
//check which of the MEs exist and respond appropriately | ||
if (!eff_v_lumi_forward) { | ||
edm::LogWarning("SiPixelPhase1EfficiencyExtras") | ||
<< "no hitefficiency_per_Lumisection_per_PXDisk_PXForward ME is available in " << effFolderName_ << std::endl; | ||
return; | ||
} | ||
if (!eff_v_lumi_barrel) { | ||
edm::LogWarning("SiPixelPhase1EfficiencyExtras") | ||
<< "no hitefficiency_per_Lumisection_per_PXLayer_PXBarrel ME is available in " << effFolderName_ << std::endl; | ||
return; | ||
} | ||
if (!vtx_v_lumi) { | ||
edm::LogWarning("SiPixelPhase1EfficiencyExtras") | ||
<< "no NumberOfGoodPVtxVsLS_GenTK ME is available in " << vtxFolderName_ << std::endl; | ||
createNvtx = false; | ||
} | ||
if (!scalLumi_v_lumi) { | ||
edm::LogWarning("SiPixelPhase1EfficiencyExtras") | ||
<< "no lumiVsLS ME is available in " << instLumiFolderName_ << std::endl; | ||
createInstLumi = false; | ||
} | ||
|
||
//If the existing MEs are empty, set the boolean to skip booking | ||
if (vtx_v_lumi && vtx_v_lumi->getEntries() == 0) | ||
createNvtx = false; | ||
if (scalLumi_v_lumi && scalLumi_v_lumi->getEntries() == 0) | ||
createInstLumi = false; | ||
|
||
double eff = 0.0; | ||
|
||
//Will pass if nvtx ME exists and is not empty | ||
if (createNvtx) { | ||
//Book new histos | ||
MonitorElement* eff_v_vtx_barrel = | ||
iBooker.book2D("hitefficiency_per_meanNvtx_per_PXLayer_PXBarrel", | ||
"hitefficiency_per_meanNvtx_per_PXLayer_PXBarrel; meanNvtx; PXLayer", | ||
500, | ||
0, | ||
100, | ||
3, | ||
.5, | ||
3.5); | ||
|
||
MonitorElement* eff_v_vtx_forward = | ||
iBooker.book2D("hitefficiency_per_meanNvtx_per_PXDisk_PXForward", | ||
"hitefficiency_per_meanNvtx_per_PXDisk_PXForward; meanNvtx; PXDisk", | ||
500, | ||
0, | ||
100, | ||
7, | ||
-3.5, | ||
3.5); | ||
|
||
//initialize variables | ||
int numLumiNvtx = int(vtx_v_lumi->getNbinsX()); | ||
double nvtx = 0.0; | ||
int binNumVtx = 0; | ||
|
||
//For loop to loop through lumisections | ||
for (int iLumi = 1; iLumi < numLumiNvtx - 1; iLumi++) { | ||
//get the meanNvtx for each lumi | ||
nvtx = vtx_v_lumi->getBinContent(iLumi); | ||
|
||
//Filter out useless iterations | ||
if (nvtx != 0) { | ||
//Grab the bin number for the nvtx | ||
binNumVtx = eff_v_vtx_barrel->getTH2F()->FindBin(nvtx); | ||
|
||
//loop through the layers | ||
for (int iLayer = 1; iLayer < 8; iLayer++) { | ||
//get the eff at the lumisection and layer | ||
eff = eff_v_lumi_forward->getBinContent(iLumi - 1, iLayer); | ||
|
||
//set the efficiency in the new histo | ||
eff_v_vtx_forward->setBinContent(binNumVtx, iLayer, eff); | ||
} | ||
|
||
//loop through the layers | ||
for (int iLayer = 1; iLayer < 5; iLayer++) { | ||
//get the efficiency for each lumi at each layer | ||
eff = eff_v_lumi_barrel->getBinContent(iLumi - 1, iLayer); | ||
|
||
//set the efficiency | ||
eff_v_vtx_barrel->setBinContent(binNumVtx, iLayer, eff); | ||
} | ||
} | ||
} | ||
} | ||
// Will pass if InstLumi ME exists and is not empty | ||
if (createInstLumi) { | ||
//Get the max value of inst lumi for plot | ||
int yMax2 = scalLumi_v_lumi->getTProfile()->GetMaximum(); | ||
yMax2 = yMax2 + yMax2 * .1; | ||
|
||
//Book new histos | ||
MonitorElement* eff_v_scalLumi_barrel = | ||
iBooker.book2D("hitefficiency_per_scalLumi_per_PXLayer_PXBarrel", | ||
"hitefficiency_per_scalLumi_per_PXLayer_PXBarrel; scal inst lumi E30; PXLayer", | ||
500, | ||
0, | ||
yMax2, | ||
3, | ||
.5, | ||
3.5); | ||
|
||
MonitorElement* eff_v_scalLumi_forward = | ||
iBooker.book2D("hitefficiency_per_scalLumi_per_PXDisk_PXForward", | ||
"hitefficiency_per_scalLumi_per_PXDisk_PXForward; scal inst lumi E30; PXDisk", | ||
500, | ||
0, | ||
yMax2, | ||
7, | ||
-3.5, | ||
3.5); | ||
|
||
//initialize variables | ||
int numLumiScal = int(scalLumi_v_lumi->getNbinsX()); | ||
double scalLumi = 0.0; | ||
int binNumScal = 0; | ||
|
||
//For loop to loop through lumisections | ||
for (int iLumi = 1; iLumi < numLumiScal - 1; iLumi++) { | ||
//get the inst lumi for each lumi | ||
scalLumi = scalLumi_v_lumi->getBinContent(iLumi); | ||
|
||
//Filter out useless iterations | ||
if (scalLumi != 0) { | ||
//Grab the bin number for the inst lumi | ||
binNumScal = eff_v_scalLumi_barrel->getTH2F()->FindBin(scalLumi); | ||
|
||
//loop through the layers | ||
for (int iLayer = 1; iLayer < 8; iLayer++) { | ||
//get the eff at the lumisection and layer | ||
eff = eff_v_lumi_forward->getBinContent(iLumi - 1, iLayer); | ||
|
||
//set the efficiency in the new histo | ||
eff_v_scalLumi_forward->setBinContent(binNumScal, iLayer, eff); | ||
} | ||
|
||
//loop through the layers | ||
for (int iLayer = 1; iLayer < 5; iLayer++) { | ||
//get the eff at the lumisection and layer | ||
eff = eff_v_lumi_barrel->getBinContent(iLumi - 1, iLayer); | ||
|
||
//set the efficiency in the new histo | ||
eff_v_scalLumi_barrel->setBinContent(binNumScal, iLayer, eff); | ||
} | ||
} | ||
} | ||
} else | ||
return; | ||
} | ||
|
||
//define this as a plug-in | ||
DEFINE_FWK_MODULE(SiPixelPhase1EfficiencyExtras); |
9 changes: 9 additions & 0 deletions
9
DQM/SiPixelPhase1Track/python/SiPixelPhase1EfficiencyExtras_cfi.py
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,9 @@ | ||
import FWCore.ParameterSet.Config as cms | ||
from DQMServices.Core.DQMEDHarvester import DQMEDHarvester | ||
|
||
|
||
SiPixelPhase1EfficiencyExtras = DQMEDHarvester("SiPixelPhase1EfficiencyExtras", | ||
EffFolderName = cms.string('PixelPhase1/Tracks/'), | ||
VtxFolderName = cms.string('Tracking/TrackParameters/generalTracks/GeneralProperties/'), | ||
InstLumiFolderName = cms.string('HLT/LumiMonitoring/') | ||
) |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters