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 
26 }
27 
28 
30 {
31  ibooker_.cd();;
32  ibooker_.setCurrentFolder("HLT/HI/" + triggerPath_);
33 
34  h_SumEt = ibooker_.book1D("h_SumEt", "SumEt", nEt, minEt, maxEt);
35  h_SiPixelClusters = ibooker_.book1D("h_SiPixelClusters", "h_SiPixelClusters", nClusters, minClusters, maxClusters);
36  h_SumEt_SiPixelClusters = ibooker_.book2D("h_SumEt_SiPixelClusters","h_SumEt_SiPixelClusters",nEt, minEt, maxEt, nClusters, minClusters, maxClusters);
37 
38  ibooker_.cd();
39 }
40 
42 {
43 
45  e.getByToken(triggerResults_,hltresults);
46  if(!hltresults.isValid())
47  {
48  edm::LogError ("HeavyIonUCCDQM") << "invalid collection: TriggerResults" << "\n";
49  return;
50  }
51 
52  bool hasFired = false;
53  const edm::TriggerNames& trigNames = e.triggerNames(*hltresults);
54  unsigned int numTriggers = trigNames.size();
55  for( unsigned int hltIndex=0; hltIndex<numTriggers; ++hltIndex ) {
56  if (trigNames.triggerName(hltIndex).find(triggerPath_) != std::string::npos && hltresults->wasrun(hltIndex) && hltresults->accept(hltIndex)){
57  hasFired = true;
58  }
59  }
60 
61  if (!hasFired) return;
62 
64  e.getByToken(theSiPixelCluster, cluster);
65  if ( cluster.isValid() ) {
66  h_SiPixelClusters->Fill(cluster->dataSize());
67  }
68 
70  e.getByToken(theCaloMet, calometColl);
71  if ( calometColl.isValid() ) {
72  h_SumEt->Fill((*calometColl).front().sumEt());
73  }
74 
75  if ( cluster.isValid() && calometColl.isValid() ) {
76  h_SumEt_SiPixelClusters->Fill((*calometColl).front().sumEt(), cluster->dataSize());
77  }
78 
79 
80 }
81 
T getParameter(std::string const &) const
bool wasrun() const
Was at least one path run?
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:234
MonitorElement * h_SiPixelClusters
MonitorElement * h_SumEt
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
void cd(void)
Definition: DQMStore.cc:269
#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:115
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override
bool isValid() const
Definition: HandleBase.h:75
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
virtual ~HeavyIonUCCDQM()
MonitorElement * h_SumEt_SiPixelClusters
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: Run.h:42