CMS 3D CMS Logo

HeavyIonUCCDQM.cc
Go to the documentation of this file.
7 
9  triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("triggerResults"));
10  theCaloMet = consumes<reco::CaloMETCollection>(ps.getParameter<edm::InputTag>("caloMet"));
11  theSiPixelCluster = consumes<edmNew::DetSetVector<SiPixelCluster> >(ps.getParameter<edm::InputTag>("pixelCluster"));
12  triggerPath_ = ps.getParameter<std::string>("triggerPath");
13 
14  nClusters = ps.getParameter<int>("nClusters");
15  minClusters = ps.getParameter<int>("minClusters");
16  maxClusters = ps.getParameter<int>("maxClusters");
17  nEt = ps.getParameter<int>("nEt");
18  minEt = ps.getParameter<double>("minEt");
19  maxEt = ps.getParameter<double>("maxEt");
20 }
21 
23 
25  ibooker_.cd();
26  ;
27  ibooker_.setCurrentFolder("HLT/HI/" + triggerPath_);
28 
29  h_SumEt = ibooker_.book1D("h_SumEt", "SumEt", nEt, minEt, maxEt);
30  h_SiPixelClusters = ibooker_.book1D("h_SiPixelClusters", "h_SiPixelClusters", nClusters, minClusters, maxClusters);
32  "h_SumEt_SiPixelClusters", "h_SumEt_SiPixelClusters", nEt, minEt, maxEt, nClusters, minClusters, maxClusters);
33 
34  ibooker_.cd();
35 }
36 
37 void HeavyIonUCCDQM::analyze(edm::Event const& e, edm::EventSetup const& eSetup) {
39  e.getByToken(triggerResults_, hltresults);
40  if (!hltresults.isValid()) {
41  edm::LogError("HeavyIonUCCDQM") << "invalid collection: TriggerResults"
42  << "\n";
43  return;
44  }
45 
46  bool hasFired = false;
47  const edm::TriggerNames& trigNames = e.triggerNames(*hltresults);
48  unsigned int numTriggers = trigNames.size();
49  for (unsigned int hltIndex = 0; hltIndex < numTriggers; ++hltIndex) {
50  if (trigNames.triggerName(hltIndex).find(triggerPath_) != std::string::npos && hltresults->wasrun(hltIndex) &&
51  hltresults->accept(hltIndex)) {
52  hasFired = true;
53  }
54  }
55 
56  if (!hasFired)
57  return;
58 
60  e.getByToken(theSiPixelCluster, cluster);
61  if (cluster.isValid()) {
62  h_SiPixelClusters->Fill(cluster->dataSize());
63  }
64 
66  e.getByToken(theCaloMet, calometColl);
67  if (calometColl.isValid()) {
68  h_SumEt->Fill((*calometColl).front().sumEt());
69  }
70 
71  if (cluster.isValid() && calometColl.isValid()) {
72  h_SumEt_SiPixelClusters->Fill((*calometColl).front().sumEt(), cluster->dataSize());
73  }
74 }
75 
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
~HeavyIonUCCDQM() override
MonitorElement * h_SiPixelClusters
MonitorElement * h_SumEt
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< edm::TriggerResults > triggerResults_
edm::EDGetTokenT< reco::CaloMETCollection > theCaloMet
std::string triggerPath_
Log< level::Error, false > LogError
void Fill(long long x)
HeavyIonUCCDQM(const edm::ParameterSet &ps)
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > theSiPixelCluster
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:57
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
bool isValid() const
Definition: HandleBase.h:70
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * h_SumEt_SiPixelClusters
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: Run.h:45