CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
DQMExample_Step1 Class Reference

#include <DQMExample_Step1.h>

Inheritance diagram for DQMExample_Step1:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

Public Member Functions

 DQMExample_Step1 (const edm::ParameterSet &ps)
 
 ~DQMExample_Step1 () override
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Protected Member Functions

void analyze (edm::Event const &e, edm::EventSetup const &eSetup) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (edm::Run const &, edm::EventSetup const &) override
 
- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 

Private Member Functions

void bookHistos (DQMStore::IBooker &)
 
double calcDeltaPhi (double phi1, double phi2)
 
double Distance (const reco::Candidate &c1, const reco::Candidate &c2)
 
double DistancePhi (const reco::Candidate &c1, const reco::Candidate &c2)
 
bool MediumEle (const edm::Event &iEvent, const edm::EventSetup &iESetup, const reco::GsfElectron &electron)
 

Private Attributes

MonitorElementh_eEta_leading
 
MonitorElementh_eEta_leading_HLT
 
MonitorElementh_eEta_leading_HLT_matched
 
MonitorElementh_eEta_leading_matched
 
MonitorElementh_eMultiplicity
 
MonitorElementh_eMultiplicity_HLT
 
MonitorElementh_ePhi_leading
 
MonitorElementh_ePhi_leading_HLT
 
MonitorElementh_ePhi_leading_HLT_matched
 
MonitorElementh_ePhi_leading_matched
 
MonitorElementh_ePt_diff
 
MonitorElementh_ePt_leading
 
MonitorElementh_ePt_leading_HLT
 
MonitorElementh_ePt_leading_HLT_matched
 
MonitorElementh_ePt_leading_matched
 
MonitorElementh_jEta_leading
 
MonitorElementh_jMultiplicity
 
MonitorElementh_jPhi_leading
 
MonitorElementh_jPt_leading
 
MonitorElementh_pfMet
 
MonitorElementh_vertex_number
 
int nBJets
 
int nElectrons
 
double ptThrJet_
 
double ptThrL1_
 
double ptThrL2_
 
double ptThrMet_
 
math::XYZPoint PVPoint_
 
edm::EDGetTokenT< reco::BeamSpottheBSCollection_
 
edm::EDGetTokenT< reco::CaloJetCollectiontheCaloJetCollection_
 
edm::EDGetTokenT< reco::ConversionCollectiontheConversionCollection_
 
edm::EDGetTokenT< reco::GsfElectronCollectiontheElectronCollection_
 
edm::EDGetTokenT< reco::PFMETCollectionthePfMETCollection_
 
edm::EDGetTokenT< reco::VertexCollectionthePVCollection_
 
edm::EDGetTokenT< trigger::TriggerEventtriggerEvent_
 
edm::InputTag triggerFilter_
 
std::string triggerPath_
 
edm::EDGetTokenT< edm::TriggerResultstriggerResults_
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr< DQMEDAnalyzerGlobalCacheinitializeGlobalCache (edm::ParameterSet const &)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

Definition at line 46 of file DQMExample_Step1.h.

Constructor & Destructor Documentation

◆ DQMExample_Step1()

DQMExample_Step1::DQMExample_Step1 ( const edm::ParameterSet ps)

Definition at line 21 of file DQMExample_Step1.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), ptThrJet_, ptThrL1_, ptThrL2_, ptThrMet_, AlCaHLTBitMon_QueryRunRegistry::string, theBSCollection_, theCaloJetCollection_, theConversionCollection_, theElectronCollection_, thePfMETCollection_, thePVCollection_, triggerEvent_, triggerFilter_, triggerPath_, and triggerResults_.

21  {
22  edm::LogInfo("DQMExample_Step1") << "Constructor DQMExample_Step1::DQMExample_Step1 " << std::endl;
23 
24  // Get parameters from configuration file
25  theElectronCollection_ = consumes<reco::GsfElectronCollection>(ps.getParameter<edm::InputTag>("electronCollection"));
26  theCaloJetCollection_ = consumes<reco::CaloJetCollection>(ps.getParameter<edm::InputTag>("caloJetCollection"));
27  thePfMETCollection_ = consumes<reco::PFMETCollection>(ps.getParameter<edm::InputTag>("pfMETCollection"));
29  consumes<reco::ConversionCollection>(ps.getParameter<edm::InputTag>("conversionsCollection"));
30  thePVCollection_ = consumes<reco::VertexCollection>(ps.getParameter<edm::InputTag>("PVCollection"));
31  theBSCollection_ = consumes<reco::BeamSpot>(ps.getParameter<edm::InputTag>("beamSpotCollection"));
32  triggerEvent_ = consumes<trigger::TriggerEvent>(ps.getParameter<edm::InputTag>("TriggerEvent"));
33  triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("TriggerResults"));
34  triggerFilter_ = ps.getParameter<edm::InputTag>("TriggerFilter");
35  triggerPath_ = ps.getParameter<std::string>("TriggerPath");
36 
37  // cuts:
38  ptThrL1_ = ps.getUntrackedParameter<double>("PtThrL1");
39  ptThrL2_ = ps.getUntrackedParameter<double>("PtThrL2");
40  ptThrJet_ = ps.getUntrackedParameter<double>("PtThrJet");
41  ptThrMet_ = ps.getUntrackedParameter<double>("PtThrMet");
42 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< reco::GsfElectronCollection > theElectronCollection_
edm::EDGetTokenT< reco::ConversionCollection > theConversionCollection_
edm::EDGetTokenT< trigger::TriggerEvent > triggerEvent_
edm::EDGetTokenT< reco::BeamSpot > theBSCollection_
T getUntrackedParameter(std::string const &, T const &) const
Log< level::Info, false > LogInfo
edm::EDGetTokenT< edm::TriggerResults > triggerResults_
edm::EDGetTokenT< reco::PFMETCollection > thePfMETCollection_
edm::InputTag triggerFilter_
std::string triggerPath_
edm::EDGetTokenT< reco::VertexCollection > thePVCollection_
edm::EDGetTokenT< reco::CaloJetCollection > theCaloJetCollection_

◆ ~DQMExample_Step1()

DQMExample_Step1::~DQMExample_Step1 ( )
override

Definition at line 47 of file DQMExample_Step1.cc.

47  {
48  edm::LogInfo("DQMExample_Step1") << "Destructor DQMExample_Step1::~DQMExample_Step1 " << std::endl;
49 }
Log< level::Info, false > LogInfo

Member Function Documentation

◆ analyze()

void DQMExample_Step1::analyze ( edm::Event const &  e,
edm::EventSetup const &  eSetup 
)
overrideprotectedvirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 72 of file DQMExample_Step1.cc.

References funct::abs(), pdwgLeptonRecoSkim_cfi::caloJetCollection, PbPb_ZMuSkimMuonDPG_cff::deltaR, Distance(), MillePedeFileConverter_cfg::e, pdwgLeptonRecoSkim_cfi::electronCollection, reco::Particle::eta(), reco::LeafCandidate::eta(), dqm::impl::MonitorElement::Fill(), h_eEta_leading, h_eEta_leading_HLT, h_eEta_leading_HLT_matched, h_eEta_leading_matched, h_eMultiplicity, h_eMultiplicity_HLT, h_ePhi_leading, h_ePhi_leading_HLT, h_ePhi_leading_HLT_matched, h_ePhi_leading_matched, h_ePt_diff, h_ePt_leading, h_ePt_leading_HLT, h_ePt_leading_HLT_matched, h_ePt_leading_matched, h_jEta_leading, h_jMultiplicity, h_jPhi_leading, h_jPt_leading, h_pfMet, h_vertex_number, HLTBitAnalyser_cfi::hltresults, edm::HandleBase::isValid(), dqmiolumiharvest::j, relativeConstraints::keys, MediumEle(), trigger::TriggerObject::particle(), reco::Particle::pdgId(), B2GDQM_cfi::pfMETCollection, reco::Particle::phi(), reco::LeafCandidate::phi(), reco::Particle::pt(), reco::LeafCandidate::pt(), ptThrJet_, ptThrL1_, ptThrL2_, PVPoint_, theCaloJetCollection_, theElectronCollection_, thePfMETCollection_, thePVCollection_, PDWG_DiPhoton_SD_cff::triggerEvent, triggerEvent_, triggerFilter_, triggerMatchMonitor_cfi::triggerObjects, triggerPath_, triggerResults_, trigNames, and findQualityFiles::v.

72  {
73  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::analyze" << std::endl;
74 
75  //-------------------------------
76  //--- Vertex Info
77  //-------------------------------
79  e.getByToken(thePVCollection_, vertexHandle);
80  if (!vertexHandle.isValid()) {
81  edm::LogError("DQMClientExample") << "invalid collection: vertex"
82  << "\n";
83  return;
84  }
85 
86  int vertex_number = vertexHandle->size();
87  reco::VertexCollection::const_iterator v = vertexHandle->begin();
88 
89  math::XYZPoint PVPoint(-999, -999, -999);
90  if (vertex_number != 0)
91  PVPoint = math::XYZPoint(v->position().x(), v->position().y(), v->position().z());
92 
93  PVPoint_ = PVPoint;
94 
95  //-------------------------------
96  //--- MET
97  //-------------------------------
100  if (!pfMETCollection.isValid()) {
101  edm::LogError("DQMClientExample") << "invalid collection: MET"
102  << "\n";
103  return;
104  }
105  //-------------------------------
106  //--- Electrons
107  //-------------------------------
110  if (!electronCollection.isValid()) {
111  edm::LogError("DQMClientExample") << "invalid collection: electrons"
112  << "\n";
113  return;
114  }
115 
116  float nEle = 0;
117  int posEle = 0, negEle = 0;
118  const reco::GsfElectron *ele1 = nullptr;
119  const reco::GsfElectron *ele2 = nullptr;
120  for (reco::GsfElectronCollection::const_iterator recoElectron = electronCollection->begin();
121  recoElectron != electronCollection->end();
122  ++recoElectron) {
123  // decreasing pT
124  if (MediumEle(e, eSetup, *recoElectron)) {
125  if (!ele1 && recoElectron->pt() > ptThrL1_)
126  ele1 = &(*recoElectron);
127 
128  else if (!ele2 && recoElectron->pt() > ptThrL2_)
129  ele2 = &(*recoElectron);
130  }
131 
132  if (recoElectron->charge() == 1)
133  posEle++;
134  else if (recoElectron->charge() == -1)
135  negEle++;
136 
137  } // end of loop over electrons
138 
139  nEle = posEle + negEle;
140 
141  //-------------------------------
142  //--- Jets
143  //-------------------------------
146  if (!caloJetCollection.isValid()) {
147  edm::LogError("DQMClientExample") << "invalid collection: jets"
148  << "\n";
149  return;
150  }
151 
152  int nJet = 0;
153  const reco::CaloJet *jet1 = nullptr;
154  const reco::CaloJet *jet2 = nullptr;
155 
156  for (reco::CaloJetCollection::const_iterator i_calojet = caloJetCollection->begin();
157  i_calojet != caloJetCollection->end();
158  ++i_calojet) {
159  // remove jet-ele matching
160  if (ele1)
161  if (Distance(*i_calojet, *ele1) < 0.3)
162  continue;
163 
164  if (ele2)
165  if (Distance(*i_calojet, *ele2) < 0.3)
166  continue;
167 
168  if (i_calojet->pt() < ptThrJet_)
169  continue;
170 
171  nJet++;
172 
173  if (!jet1)
174  jet1 = &(*i_calojet);
175 
176  else if (!jet2)
177  jet2 = &(*i_calojet);
178  }
179 
180  // ---------------------------
181  // ---- Analyze Trigger Event
182  // ---------------------------
183 
184  // check what is in the menu
186  e.getByToken(triggerResults_, hltresults);
187 
188  if (!hltresults.isValid()) {
189  edm::LogError("DQMClientExample") << "invalid collection: TriggerResults"
190  << "\n";
191  return;
192  }
193 
194  bool hasFired = false;
195  const edm::TriggerNames &trigNames = e.triggerNames(*hltresults);
196  unsigned int numTriggers = trigNames.size();
197 
198  for (unsigned int hltIndex = 0; hltIndex < numTriggers; ++hltIndex) {
199  if (trigNames.triggerName(hltIndex) == triggerPath_ && hltresults->wasrun(hltIndex) && hltresults->accept(hltIndex))
200  hasFired = true;
201  }
202 
203  // access the trigger event
205  e.getByToken(triggerEvent_, triggerEvent);
206  if (triggerEvent.failedToGet()) {
207  edm::LogError("DQMClientExample") << "invalid collection: TriggerEvent"
208  << "\n";
209  return;
210  }
211 
212  reco::Particle *ele1_HLT = nullptr;
213  int nEle_HLT = 0;
214 
215  size_t filterIndex = triggerEvent->filterIndex(triggerFilter_);
217  if (!(filterIndex >= triggerEvent->sizeFilters())) {
218  const trigger::Keys &keys = triggerEvent->filterKeys(filterIndex);
219  std::vector<reco::Particle> triggeredEle;
220 
221  for (size_t j = 0; j < keys.size(); ++j) {
223  if (abs(foundObject.particle().pdgId()) != 11)
224  continue; // make sure that it is an electron
225 
226  triggeredEle.push_back(foundObject.particle());
227  ++nEle_HLT;
228  }
229 
230  if (!triggeredEle.empty())
231  ele1_HLT = &(triggeredEle.at(0));
232  }
233 
234  //-------------------------------
235  //--- Fill the histos
236  //-------------------------------
237 
238  // vertex
239  h_vertex_number->Fill(vertex_number);
240 
241  // met
242  h_pfMet->Fill(pfMETCollection->begin()->et());
243 
244  // multiplicities
245  h_eMultiplicity->Fill(nEle);
246  h_jMultiplicity->Fill(nJet);
247  h_eMultiplicity_HLT->Fill(nEle_HLT);
248 
249  // leading not matched
250  if (ele1) {
251  h_ePt_leading->Fill(ele1->pt());
252  h_eEta_leading->Fill(ele1->eta());
253  h_ePhi_leading->Fill(ele1->phi());
254  }
255  if (ele1_HLT) {
256  h_ePt_leading_HLT->Fill(ele1_HLT->pt());
257  h_eEta_leading_HLT->Fill(ele1_HLT->eta());
258  h_ePhi_leading_HLT->Fill(ele1_HLT->phi());
259  }
260  // leading Jet
261  if (jet1) {
262  h_jPt_leading->Fill(jet1->pt());
263  h_jEta_leading->Fill(jet1->eta());
264  h_jPhi_leading->Fill(jet1->phi());
265  }
266 
267  // fill only when the trigger candidate mathes with the reco one
268  if (ele1 && ele1_HLT && deltaR(*ele1_HLT, *ele1) < 0.3 && hasFired == true) {
269  h_ePt_leading_matched->Fill(ele1->pt());
270  h_eEta_leading_matched->Fill(ele1->eta());
271  h_ePhi_leading_matched->Fill(ele1->phi());
272 
273  h_ePt_leading_HLT_matched->Fill(ele1_HLT->pt());
274  h_eEta_leading_HLT_matched->Fill(ele1_HLT->eta());
275  h_ePhi_leading_HLT_matched->Fill(ele1_HLT->phi());
276 
277  h_ePt_diff->Fill(ele1->pt() - ele1_HLT->pt());
278  }
279 }
MonitorElement * h_ePhi_leading_matched
MonitorElement * h_jMultiplicity
math::XYZPoint PVPoint_
bool MediumEle(const edm::Event &iEvent, const edm::EventSetup &iESetup, const reco::GsfElectron &electron)
edm::EDGetTokenT< reco::GsfElectronCollection > theElectronCollection_
Jets made from CaloTowers.
Definition: CaloJet.h:27
double pt() const final
transverse momentum
MonitorElement * h_ePhi_leading_HLT_matched
MonitorElement * h_ePhi_leading
MonitorElement * h_eEta_leading_HLT
MonitorElement * h_jPhi_leading
edm::EDGetTokenT< trigger::TriggerEvent > triggerEvent_
Log< level::Error, false > LogError
double Distance(const reco::Candidate &c1, const reco::Candidate &c2)
MonitorElement * h_eEta_leading
MonitorElement * h_ePt_leading_HLT
void Fill(long long x)
Single trigger physics object (e.g., an isolated muon)
Definition: TriggerObject.h:21
MonitorElement * h_ePt_diff
MonitorElement * h_ePt_leading_HLT_matched
MonitorElement * h_jEta_leading
int pdgId() const
PDG identifier.
Definition: Particle.h:150
double phi() const
momentum azimuthal angle
Definition: Particle.h:122
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorElement * h_eEta_leading_HLT_matched
MonitorElement * h_vertex_number
MonitorElement * h_ePt_leading
MonitorElement * h_eMultiplicity
MonitorElement * h_pfMet
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:75
Log< level::Info, false > LogInfo
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:57
double pt() const
transverse momentum
Definition: Particle.h:120
edm::EDGetTokenT< edm::TriggerResults > triggerResults_
MonitorElement * h_jPt_leading
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< size_type > Keys
MonitorElement * h_ePt_leading_matched
edm::EDGetTokenT< reco::PFMETCollection > thePfMETCollection_
edm::InputTag triggerFilter_
bool isValid() const
Definition: HandleBase.h:70
reco::Particle particle(reco::Particle::Charge q=0, const reco::Particle::Point &vertex=reco::Particle::Point(0, 0, 0), int status=0, bool integerCharge=true) const
Definition: TriggerObject.h:65
MonitorElement * h_ePhi_leading_HLT
std::string triggerPath_
double eta() const
momentum pseudorapidity
Definition: Particle.h:126
MonitorElement * h_eMultiplicity_HLT
double phi() const final
momentum azimuthal angle
MonitorElement * h_eEta_leading_matched
edm::EDGetTokenT< reco::VertexCollection > thePVCollection_
edm::EDGetTokenT< reco::CaloJetCollection > theCaloJetCollection_
double eta() const final
momentum pseudorapidity

◆ bookHistograms()

void DQMExample_Step1::bookHistograms ( DQMStore::IBooker ibooker_,
edm::Run const &  ,
edm::EventSetup const &   
)
overrideprotectedvirtual

Implements DQMEDAnalyzer.

Definition at line 62 of file DQMExample_Step1.cc.

References bookHistos().

62  {
63  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::bookHistograms" << std::endl;
64 
65  // book at beginRun
66  bookHistos(ibooker_);
67 }
void bookHistos(DQMStore::IBooker &)
Log< level::Info, false > LogInfo

◆ bookHistos()

void DQMExample_Step1::bookHistos ( DQMStore::IBooker ibooker_)
private

Definition at line 289 of file DQMExample_Step1.cc.

References dqm::implementation::IBooker::book1D(), dqm::implementation::NavigatorBase::cd(), h_eEta_leading, h_eEta_leading_HLT, h_eEta_leading_HLT_matched, h_eEta_leading_matched, h_eMultiplicity, h_eMultiplicity_HLT, h_ePhi_leading, h_ePhi_leading_HLT, h_ePhi_leading_HLT_matched, h_ePhi_leading_matched, h_ePt_diff, h_ePt_leading, h_ePt_leading_HLT, h_ePt_leading_HLT_matched, h_ePt_leading_matched, h_jEta_leading, h_jMultiplicity, h_jPhi_leading, h_jPt_leading, h_pfMet, h_vertex_number, and dqm::implementation::NavigatorBase::setCurrentFolder().

Referenced by bookHistograms().

289  {
290  ibooker_.cd();
291  ibooker_.setCurrentFolder("Physics/TopTest");
292 
293  h_vertex_number = ibooker_.book1D("Vertex_number", "Number of event vertices in collection", 40, -0.5, 39.5);
294 
295  h_pfMet = ibooker_.book1D("pfMet", "Pf Missing E_{T}; GeV", 20, 0.0, 100);
296 
297  h_eMultiplicity = ibooker_.book1D("NElectrons", "# of electrons per event", 10, 0., 10.);
298  h_ePt_leading_matched = ibooker_.book1D("ElePt_leading_matched", "Pt of leading electron", 50, 0., 100.);
299  h_eEta_leading_matched = ibooker_.book1D("EleEta_leading_matched", "Eta of leading electron", 50, -5., 5.);
300  h_ePhi_leading_matched = ibooker_.book1D("ElePhi_leading_matched", "Phi of leading electron", 50, -3.5, 3.5);
301 
302  h_ePt_leading = ibooker_.book1D("ElePt_leading", "Pt of leading electron", 50, 0., 100.);
303  h_eEta_leading = ibooker_.book1D("EleEta_leading", "Eta of leading electron", 50, -5., 5.);
304  h_ePhi_leading = ibooker_.book1D("ElePhi_leading", "Phi of leading electron", 50, -3.5, 3.5);
305 
306  h_jMultiplicity = ibooker_.book1D("NJets", "# of electrons per event", 10, 0., 10.);
307  h_jPt_leading = ibooker_.book1D("JetPt_leading", "Pt of leading Jet", 150, 0., 300.);
308  h_jEta_leading = ibooker_.book1D("JetEta_leading", "Eta of leading Jet", 50, -5., 5.);
309  h_jPhi_leading = ibooker_.book1D("JetPhi_leading", "Phi of leading Jet", 50, -3.5, 3.5);
310 
311  h_eMultiplicity_HLT = ibooker_.book1D("NElectrons_HLT", "# of electrons per event @HLT", 10, 0., 10.);
312  h_ePt_leading_HLT = ibooker_.book1D("ElePt_leading_HLT", "Pt of leading electron @HLT", 50, 0., 100.);
313  h_eEta_leading_HLT = ibooker_.book1D("EleEta_leading_HLT", "Eta of leading electron @HLT", 50, -5., 5.);
314  h_ePhi_leading_HLT = ibooker_.book1D("ElePhi_leading_HLT", "Phi of leading electron @HLT", 50, -3.5, 3.5);
315 
316  h_ePt_leading_HLT_matched = ibooker_.book1D("ElePt_leading_HLT_matched", "Pt of leading electron @HLT", 50, 0., 100.);
318  ibooker_.book1D("EleEta_leading_HLT_matched", "Eta of leading electron @HLT", 50, -5., 5.);
320  ibooker_.book1D("ElePhi_leading_HLT_matched", "Phi of leading electron @HLT", 50, -3.5, 3.5);
321 
322  h_ePt_diff = ibooker_.book1D("ElePt_diff_matched", "pT(RECO) - pT(HLT) for mathed candidates", 100, -10, 10.);
323 
324  ibooker_.cd();
325 }
MonitorElement * h_ePhi_leading_matched
MonitorElement * h_jMultiplicity
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
MonitorElement * h_ePhi_leading_HLT_matched
MonitorElement * h_ePhi_leading
MonitorElement * h_eEta_leading_HLT
MonitorElement * h_jPhi_leading
MonitorElement * h_eEta_leading
MonitorElement * h_ePt_leading_HLT
MonitorElement * h_ePt_diff
MonitorElement * h_ePt_leading_HLT_matched
MonitorElement * h_jEta_leading
MonitorElement * h_eEta_leading_HLT_matched
MonitorElement * h_vertex_number
MonitorElement * h_ePt_leading
MonitorElement * h_eMultiplicity
MonitorElement * h_pfMet
MonitorElement * h_jPt_leading
MonitorElement * h_ePt_leading_matched
MonitorElement * h_ePhi_leading_HLT
MonitorElement * h_eMultiplicity_HLT
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_eEta_leading_matched

◆ calcDeltaPhi()

double DQMExample_Step1::calcDeltaPhi ( double  phi1,
double  phi2 
)
private

Definition at line 338 of file DQMExample_Step1.cc.

References SiPixelRawToDigiRegional_cfi::deltaPhi.

338  {
339  double deltaPhi = phi1 - phi2;
340  if (deltaPhi < 0)
341  deltaPhi = -deltaPhi;
342  if (deltaPhi > 3.1415926) {
343  deltaPhi = 2 * 3.1415926 - deltaPhi;
344  }
345  return deltaPhi;
346 }

◆ Distance()

double DQMExample_Step1::Distance ( const reco::Candidate c1,
const reco::Candidate c2 
)
private

Definition at line 331 of file DQMExample_Step1.cc.

References alignmentValidation::c1, and PbPb_ZMuSkimMuonDPG_cff::deltaR.

Referenced by analyze().

◆ DistancePhi()

double DQMExample_Step1::DistancePhi ( const reco::Candidate c1,
const reco::Candidate c2 
)
private

Definition at line 333 of file DQMExample_Step1.cc.

References alignmentValidation::c1, SiPixelRawToDigiRegional_cfi::deltaPhi, and reco::Candidate::p4().

333  {
334  return deltaPhi(c1.p4().phi(), c2.p4().phi());
335 }
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector

◆ dqmBeginRun()

void DQMExample_Step1::dqmBeginRun ( edm::Run const &  ,
edm::EventSetup const &   
)
overrideprotectedvirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 55 of file DQMExample_Step1.cc.

55  {
56  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::beginRun" << std::endl;
57 }
Log< level::Info, false > LogInfo

◆ MediumEle()

bool DQMExample_Step1::MediumEle ( const edm::Event iEvent,
const edm::EventSetup iESetup,
const reco::GsfElectron electron 
)
private

Definition at line 352 of file DQMExample_Step1.cc.

References PVValHelper::dxy, PVValHelper::dz, HPSPFTauProducerPuppi_cfi::electron, PVValHelper::eta, ConversionTools::hasMatchedConversion(), iEvent, reco::HitPattern::MISSING_INNER_HITS, reco::BeamSpot::position(), DiDispStaMuonMonitor_cfi::pt, PVPoint_, theBSCollection_, and theConversionCollection_.

Referenced by analyze().

354  {
355  //********* CONVERSION TOOLS
357  iEvent.getByToken(theConversionCollection_, conversions_h);
358 
359  bool isMediumEle = false;
360 
361  float pt = electron.pt();
362  float eta = electron.eta();
363 
364  int isEB = electron.isEB();
365  float sigmaIetaIeta = electron.sigmaIetaIeta();
366  float DetaIn = electron.deltaEtaSuperClusterTrackAtVtx();
367  float DphiIn = electron.deltaPhiSuperClusterTrackAtVtx();
368  float HOverE = electron.hadronicOverEm();
369  float ooemoop = (1.0 / electron.ecalEnergy() - electron.eSuperClusterOverP() / electron.ecalEnergy());
370 
371  int mishits = electron.gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
372  int nAmbiguousGsfTracks = electron.ambiguousGsfTracksSize();
373 
374  reco::GsfTrackRef eleTrack = electron.gsfTrack();
375  float dxy = eleTrack->dxy(PVPoint_);
376  float dz = eleTrack->dz(PVPoint_);
377 
379  iEvent.getByToken(theBSCollection_, BSHandle);
380  const reco::BeamSpot BS = *BSHandle;
381 
382  bool isConverted = ConversionTools::hasMatchedConversion(electron, *conversions_h, BS.position());
383 
384  // default
385  if ((pt > 12.) && (fabs(eta) < 2.5) &&
386  (((isEB == 1) && (fabs(DetaIn) < 0.004)) || ((isEB == 0) && (fabs(DetaIn) < 0.007))) &&
387  (((isEB == 1) && (fabs(DphiIn) < 0.060)) || ((isEB == 0) && (fabs(DphiIn) < 0.030))) &&
388  (((isEB == 1) && (sigmaIetaIeta < 0.010)) || ((isEB == 0) && (sigmaIetaIeta < 0.030))) &&
389  (((isEB == 1) && (HOverE < 0.120)) || ((isEB == 0) && (HOverE < 0.100))) &&
390  (((isEB == 1) && (fabs(ooemoop) < 0.050)) || ((isEB == 0) && (fabs(ooemoop) < 0.050))) &&
391  (((isEB == 1) && (fabs(dxy) < 0.020)) || ((isEB == 0) && (fabs(dxy) < 0.020))) &&
392  (((isEB == 1) && (fabs(dz) < 0.100)) || ((isEB == 0) && (fabs(dz) < 0.100))) &&
393  (((isEB == 1) && (!isConverted)) || ((isEB == 0) && (!isConverted))) && (mishits == 0) &&
394  (nAmbiguousGsfTracks == 0))
395  isMediumEle = true;
396 
397  return isMediumEle;
398 }
math::XYZPoint PVPoint_
const Point & position() const
position
Definition: BeamSpot.h:59
edm::EDGetTokenT< reco::ConversionCollection > theConversionCollection_
edm::EDGetTokenT< reco::BeamSpot > theBSCollection_
int iEvent
Definition: GenABIO.cc:224
static bool hasMatchedConversion(const reco::GsfElectron &ele, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)

Member Data Documentation

◆ h_eEta_leading

MonitorElement* DQMExample_Step1::h_eEta_leading
private

Definition at line 96 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_eEta_leading_HLT

MonitorElement* DQMExample_Step1::h_eEta_leading_HLT
private

Definition at line 104 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_eEta_leading_HLT_matched

MonitorElement* DQMExample_Step1::h_eEta_leading_HLT_matched
private

Definition at line 107 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_eEta_leading_matched

MonitorElement* DQMExample_Step1::h_eEta_leading_matched
private

Definition at line 99 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_eMultiplicity

MonitorElement* DQMExample_Step1::h_eMultiplicity
private

Definition at line 94 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_eMultiplicity_HLT

MonitorElement* DQMExample_Step1::h_eMultiplicity_HLT
private

Definition at line 102 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_ePhi_leading

MonitorElement* DQMExample_Step1::h_ePhi_leading
private

Definition at line 97 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_ePhi_leading_HLT

MonitorElement* DQMExample_Step1::h_ePhi_leading_HLT
private

Definition at line 105 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_ePhi_leading_HLT_matched

MonitorElement* DQMExample_Step1::h_ePhi_leading_HLT_matched
private

Definition at line 108 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_ePhi_leading_matched

MonitorElement* DQMExample_Step1::h_ePhi_leading_matched
private

Definition at line 100 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_ePt_diff

MonitorElement* DQMExample_Step1::h_ePt_diff
private

Definition at line 115 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_ePt_leading

MonitorElement* DQMExample_Step1::h_ePt_leading
private

Definition at line 95 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_ePt_leading_HLT

MonitorElement* DQMExample_Step1::h_ePt_leading_HLT
private

Definition at line 103 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_ePt_leading_HLT_matched

MonitorElement* DQMExample_Step1::h_ePt_leading_HLT_matched
private

Definition at line 106 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_ePt_leading_matched

MonitorElement* DQMExample_Step1::h_ePt_leading_matched
private

Definition at line 98 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_jEta_leading

MonitorElement* DQMExample_Step1::h_jEta_leading
private

Definition at line 112 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_jMultiplicity

MonitorElement* DQMExample_Step1::h_jMultiplicity
private

Definition at line 110 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_jPhi_leading

MonitorElement* DQMExample_Step1::h_jPhi_leading
private

Definition at line 113 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_jPt_leading

MonitorElement* DQMExample_Step1::h_jPt_leading
private

Definition at line 111 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_pfMet

MonitorElement* DQMExample_Step1::h_pfMet
private

Definition at line 92 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ h_vertex_number

MonitorElement* DQMExample_Step1::h_vertex_number
private

Definition at line 90 of file DQMExample_Step1.h.

Referenced by analyze(), and bookHistos().

◆ nBJets

int DQMExample_Step1::nBJets
private

Definition at line 87 of file DQMExample_Step1.h.

◆ nElectrons

int DQMExample_Step1::nElectrons
private

Definition at line 86 of file DQMExample_Step1.h.

◆ ptThrJet_

double DQMExample_Step1::ptThrJet_
private

Definition at line 83 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ ptThrL1_

double DQMExample_Step1::ptThrL1_
private

Definition at line 81 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ ptThrL2_

double DQMExample_Step1::ptThrL2_
private

Definition at line 82 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ ptThrMet_

double DQMExample_Step1::ptThrMet_
private

Definition at line 84 of file DQMExample_Step1.h.

Referenced by DQMExample_Step1().

◆ PVPoint_

math::XYZPoint DQMExample_Step1::PVPoint_
private

Definition at line 67 of file DQMExample_Step1.h.

Referenced by analyze(), and MediumEle().

◆ theBSCollection_

edm::EDGetTokenT<reco::BeamSpot> DQMExample_Step1::theBSCollection_
private

Definition at line 75 of file DQMExample_Step1.h.

Referenced by DQMExample_Step1(), and MediumEle().

◆ theCaloJetCollection_

edm::EDGetTokenT<reco::CaloJetCollection> DQMExample_Step1::theCaloJetCollection_
private

Definition at line 72 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ theConversionCollection_

edm::EDGetTokenT<reco::ConversionCollection> DQMExample_Step1::theConversionCollection_
private

Definition at line 71 of file DQMExample_Step1.h.

Referenced by DQMExample_Step1(), and MediumEle().

◆ theElectronCollection_

edm::EDGetTokenT<reco::GsfElectronCollection> DQMExample_Step1::theElectronCollection_
private

Definition at line 70 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ thePfMETCollection_

edm::EDGetTokenT<reco::PFMETCollection> DQMExample_Step1::thePfMETCollection_
private

Definition at line 73 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ thePVCollection_

edm::EDGetTokenT<reco::VertexCollection> DQMExample_Step1::thePVCollection_
private

Definition at line 74 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ triggerEvent_

edm::EDGetTokenT<trigger::TriggerEvent> DQMExample_Step1::triggerEvent_
private

Definition at line 76 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ triggerFilter_

edm::InputTag DQMExample_Step1::triggerFilter_
private

Definition at line 78 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ triggerPath_

std::string DQMExample_Step1::triggerPath_
private

Definition at line 79 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().

◆ triggerResults_

edm::EDGetTokenT<edm::TriggerResults> DQMExample_Step1::triggerResults_
private

Definition at line 77 of file DQMExample_Step1.h.

Referenced by analyze(), and DQMExample_Step1().