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 */
16 #include "TPRegexp.h"
17 
18 // -----------------------------
19 // constructors and destructor
20 // -----------------------------
21 
26  consumes<SiPixelClusterCollectionNew>(parameters_.getParameter<edm::InputTag>("PixelClusterInputTag"));
28  consumes<reco::VertexCollection>(parameters_.getParameter<edm::InputTag>("PixelVertexInputTag"));
29  hltInputTagToken_ = consumes<edm::TriggerResults>(parameters_.getParameter<edm::InputTag>("HLTInputTag"));
33  minVtxDoF_ = parameters_.getParameter<double>("MinVtxDoF");
34 }
35 
37 
39  std::vector<std::string> hltPathsOfInterest =
40  parameters_.getParameter<std::vector<std::string> >("HLTPathsOfInterest");
41  if (hltPathsOfInterest.empty())
42  return;
43 
44  const std::vector<std::string>& pathList = hltConfig_.triggerNames();
45  std::vector<std::string> selectedPaths;
46  for (std::vector<std::string>::const_iterator it = pathList.begin(); it != pathList.end(); ++it) {
47  int nmatch = 0;
48  for (std::vector<std::string>::const_iterator kt = hltPathsOfInterest.begin(); kt != hltPathsOfInterest.end();
49  ++kt) {
50  nmatch += TPRegexp(*kt).Match(*it);
51  }
52  if (!nmatch)
53  continue;
54  else
55  selectedPaths.push_back(*it);
56  }
57 
58  edm::ParameterSet ClusHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1ClusPar");
59  edm::ParameterSet VtxHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1VtxPar");
60 
61  std::string currentFolder = moduleName_ + "/" + folderName_;
62  dbe_->setCurrentFolder(currentFolder);
63 
64  PixelMEs local_MEs;
65  for (std::vector<std::string>::iterator it = selectedPaths.begin(); it != selectedPaths.end(); it++) {
66  std::string tag = (*it);
67  std::map<std::string, PixelMEs>::iterator iPos = histoMap_.find(tag);
68  if (iPos == histoMap_.end()) {
69  std::string hname, htitle;
70 
71  hname = "nPxlClus_";
72  hname += tag;
73  htitle = "# of Pixel Clusters (";
74  htitle += tag + ")";
75  local_MEs.clusME = dbe_->book1D(hname,
76  htitle,
77  ClusHistoPar.getParameter<int32_t>("Xbins"),
78  ClusHistoPar.getParameter<double>("Xmin"),
79  ClusHistoPar.getParameter<double>("Xmax"));
80 
81  hname = "nPxlVtx_";
82  hname += tag;
83  htitle = "# of Pixel Vertices (";
84  htitle += tag + ")";
85  local_MEs.vtxME = dbe_->book1D(hname,
86  htitle,
87  VtxHistoPar.getParameter<int32_t>("Xbins"),
88  VtxHistoPar.getParameter<double>("Xmin"),
89  VtxHistoPar.getParameter<double>("Xmax"));
90 
91  histoMap_.insert(std::make_pair(tag, local_MEs));
92  }
93  }
94 }
95 
97 
98 void PixelVTXMonitor::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
99  bool changed = true;
100  if (hltConfig_.init(iRun, iSetup, hltInputTag_.process(), changed)) {
101  // if init returns TRUE, initialisation has succeeded!
102  edm::LogInfo("PixelVTXMonitor") << "HLT config with process name " << hltInputTag_.process()
103  << " successfully extracted";
104  } else {
105  // if init returns FALSE, initialisation has NOT succeeded, which indicates a problem
106  // with the file and/or code and needs to be investigated!
107  edm::LogError("PixelVTXMonotor") << "Error! HLT config extraction with process name " << hltInputTag_.process()
108  << " failed";
109  // In this case, all access methods will return empty values!
110  }
111  bookHistograms();
112 }
114  if (histoMap_.empty())
115  return;
116 
117  //Access Pixel Clusters
119  iEvent.getByToken(pixelClusterInputTagToken_, siPixelClusters);
120 
121  if (!siPixelClusters.isValid()) {
122  edm::LogError("PixelVTXMonotor") << "Could not find Cluster Collection " << pixelClusterInputTag_;
123  return;
124  }
125  unsigned nClusters = siPixelClusters->size();
126 
127  //Access Pixel Verteces
129  iEvent.getByToken(pixelVertexInputTagToken_, pixelVertices);
130  if (!pixelVertices.isValid()) {
131  edm::LogError("PixelVTXMonotor") << "Could not find Vertex Collection " << pixelVertexInputTag_;
132  return;
133  }
134 
135  int nVtx = 0;
136  for (reco::VertexCollection::const_iterator ivtx = pixelVertices->begin(); ivtx != pixelVertices->end(); ++ivtx) {
137  if (minVtxDoF_ == -1)
138  nVtx++;
139  else {
140  if ((ivtx->isValid() == true) && (ivtx->isFake() == false) && (ivtx->ndof() >= minVtxDoF_) &&
141  (ivtx->tracksSize() != 0))
142  nVtx++;
143  }
144  }
145  // Access Trigger Results
147  iEvent.getByToken(hltInputTagToken_, triggerResults);
148  if (!triggerResults.isValid())
149  return;
150 
151  for (std::map<std::string, PixelMEs>::iterator it = histoMap_.begin(); it != histoMap_.end(); ++it) {
152  std::string path = it->first;
153  MonitorElement* me_clus = it->second.clusME;
154  MonitorElement* me_vtx = it->second.vtxME;
155  unsigned int index = hltConfig_.triggerIndex(path);
156  if (index < triggerResults->size() && triggerResults->accept(index)) {
157  if (me_vtx)
158  me_vtx->Fill(nVtx);
159  if (me_clus)
160  me_clus->Fill(nClusters);
161  }
162  }
163 }
164 
165 void PixelVTXMonitor::endRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {}
166 
168 // Define this as a plug-in
171 
172 // Local Variables:
173 // show-trailing-whitespace: t
174 // truncate-lines: t
175 // End:
size
Write out results.
edm::EDGetTokenT< reco::VertexCollection > pixelVertexInputTagToken_
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::TriggerResults > hltInputTagToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::ParameterSet parameters_
void beginJob() override
bool accept() const
Has at least one path accepted the event?
const std::vector< std::string > & triggerNames() const
names of trigger paths
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
void endJob() override
unsigned int triggerIndex(const std::string &triggerName) const
slot position of trigger path in trigger table (0 to size-1)
void Fill(long long x)
HLTConfigProvider hltConfig_
edm::InputTag pixelClusterInputTag_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::map< std::string, PixelMEs > histoMap_
void endRun(edm::Run const &iRun, edm::EventSetup const &iSetup) override
static std::string const triggerResults
Definition: EdmProvDump.cc:45
bool isValid() const
Definition: HandleBase.h:70
edm::InputTag pixelVertexInputTag_
~PixelVTXMonitor() override
std::string folderName_
edm::InputTag hltInputTag_
edm::EDGetTokenT< SiPixelClusterCollectionNew > pixelClusterInputTagToken_
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
std::string moduleName_
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:639
std::string const & process() const
Definition: InputTag.h:40
size_type size() const
void beginRun(edm::Run const &iRun, edm::EventSetup const &iSetup) override
PixelVTXMonitor(const edm::ParameterSet &)
MonitorElement * book1D(char_string const &name, char_string const &title, int const nchX, double const lowX, double const highX)
Book 1D histogram.
Definition: DQMStore.cc:1121
Definition: Run.h:45