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