CMS 3D CMS Logo

PixelVTXMonitor.cc
Go to the documentation of this file.
1 /*
2  * \file PixelVTXMonitor.cc
3  * \author S. Dutta
4  * Last Update:
5  *
6  * Description: Pixel Vertex Monitoring for different HLT paths
7  *
8 */
9 
10 // system includes
11 #include <string>
12 #include <vector>
13 #include <map>
14 
15 // user includes
30 
31 // ROOT includes
32 #include "TPRegexp.h"
33 
34 //
35 // class declaration
36 //
37 
39 public:
43  ~PixelVTXMonitor() override = default;
44 
45 protected:
46  void bookHistograms(DQMStore::IBooker& iBooker, const edm::Run& iRun, const edm::EventSetup& iSetup) override;
47 
48 private:
49  void dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) override;
50  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
51 
53 
62  const float minVtxDoF_;
63 
65 
66  struct PixelMEs {
69  };
70 
71  std::map<std::string, PixelMEs> histoMap_;
72 };
73 
74 // -----------------------------
75 // constructors and destructor
76 // -----------------------------
77 
79  : parameters_(ps),
80  moduleName_(parameters_.getParameter<std::string>("ModuleName")),
81  folderName_(parameters_.getParameter<std::string>("FolderName")),
82  pixelClusterInputTag_(parameters_.getParameter<edm::InputTag>("PixelClusterInputTag")),
83  pixelVertexInputTag_(parameters_.getParameter<edm::InputTag>("PixelVertexInputTag")),
84  hltInputTag_(parameters_.getParameter<edm::InputTag>("HLTInputTag")),
85  pixelClusterInputTagToken_(consumes<SiPixelClusterCollectionNew>(pixelClusterInputTag_)),
86  pixelVertexInputTagToken_(consumes<reco::VertexCollection>(pixelVertexInputTag_)),
87  hltInputTagToken_(consumes<edm::TriggerResults>(hltInputTag_)),
88  minVtxDoF_(parameters_.getParameter<double>("MinVtxDoF")) {}
89 
91  std::vector<std::string> hltPathsOfInterest =
92  parameters_.getParameter<std::vector<std::string> >("HLTPathsOfInterest");
93  if (hltPathsOfInterest.empty())
94  return;
95 
96  const std::vector<std::string>& pathList = hltConfig_.triggerNames();
97  std::vector<std::string> selectedPaths;
98  for (const auto& it : pathList) {
99  int nmatch = 0;
100  for (const auto& kt : hltPathsOfInterest) {
101  nmatch += TPRegexp(kt).Match(it);
102  }
103  if (!nmatch)
104  continue;
105  else
106  selectedPaths.push_back(it);
107  }
108 
109  edm::ParameterSet ClusHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1ClusPar");
110  edm::ParameterSet VtxHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1VtxPar");
111 
112  std::string currentFolder = moduleName_ + "/" + folderName_;
113  iBooker.setCurrentFolder(currentFolder);
114 
115  PixelMEs local_MEs;
116  for (const auto& tag : selectedPaths) {
117  std::map<std::string, PixelMEs>::iterator iPos = histoMap_.find(tag);
118  if (iPos == histoMap_.end()) {
119  std::string hname, htitle;
120 
121  hname = "nPxlClus_";
122  hname += tag;
123  htitle = "# of Pixel Clusters (";
124  htitle += tag + ")";
125  local_MEs.clusME = iBooker.book1D(hname,
126  htitle,
127  ClusHistoPar.getParameter<int32_t>("Xbins"),
128  ClusHistoPar.getParameter<double>("Xmin"),
129  ClusHistoPar.getParameter<double>("Xmax"));
130 
131  hname = "nPxlVtx_";
132  hname += tag;
133  htitle = "# of Pixel Vertices (";
134  htitle += tag + ")";
135  local_MEs.vtxME = iBooker.book1D(hname,
136  htitle,
137  VtxHistoPar.getParameter<int32_t>("Xbins"),
138  VtxHistoPar.getParameter<double>("Xmin"),
139  VtxHistoPar.getParameter<double>("Xmax"));
140 
141  histoMap_.insert(std::make_pair(tag, local_MEs));
142  }
143  }
144 }
145 
146 void PixelVTXMonitor::dqmBeginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
147  bool changed = true;
148  if (hltConfig_.init(iRun, iSetup, hltInputTag_.process(), changed)) {
149  // if init returns TRUE, initialisation has succeeded!
150  edm::LogInfo("PixelVTXMonitor") << "HLT config with process name " << hltInputTag_.process()
151  << " successfully extracted";
152  } else {
153  // if init returns FALSE, initialisation has NOT succeeded, which indicates a problem
154  // with the file and/or code and needs to be investigated!
155  edm::LogError("PixelVTXMonotor") << "Error! HLT config extraction with process name " << hltInputTag_.process()
156  << " failed";
157  // In this case, all access methods will return empty values!
158  }
159 }
161  if (histoMap_.empty())
162  return;
163 
164  //Access Pixel Clusters
166 
167  if (!siPixelClusters.isValid()) {
168  edm::LogError("PixelVTXMonotor") << "Could not find Cluster Collection " << pixelClusterInputTag_;
169  return;
170  }
171  unsigned nClusters = siPixelClusters->size();
172 
173  //Access Pixel Verteces
175  if (!pixelVertices.isValid()) {
176  edm::LogError("PixelVTXMonotor") << "Could not find Vertex Collection " << pixelVertexInputTag_;
177  return;
178  }
179 
180  int nVtx = 0;
181  for (const auto& ivtx : *pixelVertices) {
182  if (minVtxDoF_ == -1)
183  nVtx++;
184  else {
185  if ((ivtx.isValid() == true) && (ivtx.isFake() == false) && (ivtx.ndof() >= minVtxDoF_) &&
186  (ivtx.tracksSize() != 0))
187  nVtx++;
188  }
189  }
190  // Access Trigger Results
192  if (!triggerResults.isValid())
193  return;
194 
195  for (const auto& it : histoMap_) {
196  std::string path = it.first;
197  MonitorElement* me_clus = it.second.clusME;
198  MonitorElement* me_vtx = it.second.vtxME;
199  unsigned int index = hltConfig_.triggerIndex(path);
200  if (index < triggerResults->size() && triggerResults->accept(index)) {
201  if (me_vtx)
202  me_vtx->Fill(nVtx);
203  if (me_clus)
204  me_clus->Fill(nClusters);
205  }
206  }
207 }
208 
209 // Define this as a plug-in
size
Write out results.
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
dqm::legacy::DQMStore DQMStore
const float minVtxDoF_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::EDGetTokenT< edm::TriggerResults > hltInputTagToken_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
dqm::legacy::MonitorElement MonitorElement
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ParameterSet parameters_
const edm::EDGetTokenT< SiPixelClusterCollectionNew > pixelClusterInputTagToken_
Log< level::Error, false > LogError
const edm::InputTag hltInputTag_
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
const edm::EDGetTokenT< reco::VertexCollection > pixelVertexInputTagToken_
void dqmBeginRun(const edm::Run &iRun, const edm::EventSetup &iSetup) override
void bookHistograms(DQMStore::IBooker &iBooker, const edm::Run &iRun, const edm::EventSetup &iSetup) override
void Fill(long long x)
HLTConfigProvider hltConfig_
int iEvent
Definition: GenABIO.cc:224
std::map< std::string, PixelMEs > histoMap_
unsigned int triggerIndex(const std::string &triggerName) const
slot position of trigger path in trigger table (0 to size-1)
static std::string const triggerResults
Definition: EdmProvDump.cc:44
const std::string folderName_
const edm::InputTag pixelClusterInputTag_
Log< level::Info, false > LogInfo
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
const std::vector< std::string > & triggerNames() const
names of trigger paths
const std::string moduleName_
PixelVTXMonitor(const edm::ParameterSet &)
fixed size matrix
HLT enums.
~PixelVTXMonitor() override=default
const edm::InputTag pixelVertexInputTag_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
std::string const & process() const
Definition: InputTag.h:40
Definition: Run.h:45