CMS 3D CMS Logo

HeavyIonUCCDQM.cc
Go to the documentation of this file.
7 
9 {
10  triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("triggerResults"));
11  theCaloMet = consumes<reco::CaloMETCollection>(ps.getParameter<edm::InputTag>("caloMet"));
12  theSiPixelCluster = consumes<edmNew::DetSetVector<SiPixelCluster> >(ps.getParameter<edm::InputTag>("pixelCluster"));
13  triggerPath_ = ps.getParameter<std::string>("triggerPath");
14 
15  nClusters = ps.getParameter<int>("nClusters");
16  minClusters = ps.getParameter<int>("minClusters");
17  maxClusters = ps.getParameter<int>("maxClusters");
18  nEt = ps.getParameter<int>("nEt");
19  minEt = ps.getParameter<double>("minEt");
20  maxEt = ps.getParameter<double>("maxEt");
21 }
22 
24 
25 
27 {
28  ibooker_.cd();;
29  ibooker_.setCurrentFolder("HLT/HI/" + triggerPath_);
30 
31  h_SumEt = ibooker_.book1D("h_SumEt", "SumEt", nEt, minEt, maxEt);
32  h_SiPixelClusters = ibooker_.book1D("h_SiPixelClusters", "h_SiPixelClusters", nClusters, minClusters, maxClusters);
33  h_SumEt_SiPixelClusters = ibooker_.book2D("h_SumEt_SiPixelClusters","h_SumEt_SiPixelClusters",nEt, minEt, maxEt, nClusters, minClusters, maxClusters);
34 
35  ibooker_.cd();
36 }
37 
39 {
40 
42  e.getByToken(triggerResults_,hltresults);
43  if(!hltresults.isValid())
44  {
45  edm::LogError ("HeavyIonUCCDQM") << "invalid collection: TriggerResults" << "\n";
46  return;
47  }
48 
49  bool hasFired = false;
50  const edm::TriggerNames& trigNames = e.triggerNames(*hltresults);
51  unsigned int numTriggers = trigNames.size();
52  for( unsigned int hltIndex=0; hltIndex<numTriggers; ++hltIndex ) {
53  if (trigNames.triggerName(hltIndex).find(triggerPath_) != std::string::npos && hltresults->wasrun(hltIndex) && hltresults->accept(hltIndex)){
54  hasFired = true;
55  }
56  }
57 
58  if (!hasFired) return;
59 
61  e.getByToken(theSiPixelCluster, cluster);
62  if ( cluster.isValid() ) {
63  h_SiPixelClusters->Fill(cluster->dataSize());
64  }
65 
67  e.getByToken(theCaloMet, calometColl);
68  if ( calometColl.isValid() ) {
69  h_SumEt->Fill((*calometColl).front().sumEt());
70  }
71 
72  if ( cluster.isValid() && calometColl.isValid() ) {
73  h_SumEt_SiPixelClusters->Fill((*calometColl).front().sumEt(), cluster->dataSize());
74  }
75 
76 
77 }
78 
T getParameter(std::string const &) const
bool wasrun() const
Was at least one path run?
~HeavyIonUCCDQM() override
MonitorElement * h_SiPixelClusters
MonitorElement * h_SumEt
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edm::TriggerResults > triggerResults_
bool accept() const
Has at least one path accepted the event?
edm::EDGetTokenT< reco::CaloMETCollection > theCaloMet
std::string triggerPath_
Strings::size_type size() const
Definition: TriggerNames.cc:39
void Fill(long long x)
HeavyIonUCCDQM(const edm::ParameterSet &ps)
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > theSiPixelCluster
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override
bool isValid() const
Definition: HandleBase.h:74
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
MonitorElement * h_SumEt_SiPixelClusters
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:301
Definition: Run.h:43