CMS 3D CMS Logo

DQMHcalIsolatedBunchAlCaReco.cc
Go to the documentation of this file.
1 /*
2  * \file DQMHcalIsolatedBunchAlCaReco.cc
3  *
4  * \author Olga Kodolova
5  *
6  *
7  *
8  * Description: Monitoring of Phi Symmetry Calibration Stream
9  */
10 
16 
19 
20 // DQM include files
21 
23 
24 // work on collections
25 
30 
36 
38 
39 // ******************************************
40 // constructors
41 // *****************************************
42 
44  //
45  // Input from configurator file
46  //
47  folderName_ = ps.getUntrackedParameter<std::string>("FolderName", "ALCAStreamHcalIsolatedBunch");
48  trigName_ = ps.getParameter<std::string>("TriggerName");
49  plotAll_ = ps.getUntrackedParameter<bool>("PlotAll", true);
50 
51  hbhereco_ = consumes<HBHERecHitCollection>(ps.getParameter<edm::InputTag>("hbheInput"));
52  horeco_ = consumes<HORecHitCollection>(ps.getParameter<edm::InputTag>("hoInput"));
53  hfreco_ = consumes<HFRecHitCollection>(ps.getParameter<edm::InputTag>("hfInput"));
54  trigResult_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("TriggerResult"));
55 }
56 
58 
59 //--------------------------------------------------------
61  edm::Run const &irun,
62  edm::EventSetup const &isetup) {
63  // create and cd into new folder
65  h_Event_ = ibooker.book1D("hEvent", "Selection summary", 10, 0, 10);
66  h_hbhehit_ = ibooker.book1D("hHBHEHit", "Size of HBHE Collection", 200, 0, 2000);
67  h_hohit_ = ibooker.book1D("hHOHit", "Size of HO Collection", 200, 0, 2000);
68  h_hfhit_ = ibooker.book1D("hHFHit", "Size of HF Collection", 200, 0, 2000);
69 }
70 
71 //-------------------------------------------------------------
72 
74  bool accept(false);
77  iEvent.getByToken(trigResult_, triggerResults);
78  if (triggerResults.isValid()) {
79  const edm::TriggerNames &triggerNames = iEvent.triggerNames(*triggerResults);
80  const std::vector<std::string> &triggerNames_ = triggerNames.triggerNames();
81  for (unsigned int iHLT = 0; iHLT < triggerResults->size(); iHLT++) {
82  int hlt = triggerResults->accept(iHLT);
83  if (hlt > 0) {
84  if (triggerNames_[iHLT].find(trigName_) != std::string::npos) {
85  accept = true;
86  break;
87  }
88  }
89  }
90  }
91  h_Event_->Fill(0.);
92  if (accept)
93  h_Event_->Fill(1.);
94 
95  if (accept || plotAll_) {
97  iEvent.getByToken(hbhereco_, hbhe);
98  if (!hbhe.isValid()) {
99  edm::LogInfo("HcalCalib") << "Cannot get hbhe product!" << std::endl;
100  } else {
101  h_hbhehit_->Fill(hbhe->size());
102  }
103 
105  iEvent.getByToken(hfreco_, hf);
106  if (!hf.isValid()) {
107  edm::LogInfo("HcalCalib") << "Cannot get hf product!" << std::endl;
108  } else {
109  h_hfhit_->Fill(hf->size());
110  }
111 
113  iEvent.getByToken(horeco_, ho);
114  if (!ho.isValid()) {
115  edm::LogInfo("HcalCalib") << "Cannot get ho product!" << std::endl;
116  } else {
117  h_hohit_->Fill(ho->size());
118  }
119  }
120 
121 } // analyze
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< HBHERecHitCollection > hbhereco_
object to monitor
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
bool accept() const
Has at least one path accepted the event?
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
Strings const & triggerNames() const
Definition: TriggerNames.cc:20
void Fill(long long x)
DQMHcalIsolatedBunchAlCaReco(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
unsigned int size() const
Get number of paths stored.
std::string folderName_
DQM folder name.
edm::EDGetTokenT< HORecHitCollection > horeco_
edm::EDGetTokenT< HFRecHitCollection > hfreco_
void analyze(const edm::Event &e, const edm::EventSetup &c) override
static std::string const triggerResults
Definition: EdmProvDump.cc:45
bool isValid() const
Definition: HandleBase.h:70
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
edm::EDGetTokenT< edm::TriggerResults > trigResult_
size_type size() const
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:265
Definition: Run.h:45