CMS 3D CMS Logo

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

#include <DrellYanValidation.h>

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

Public Member Functions

void analyze (edm::Event const &, edm::EventSetup const &) override
 
void bookHistograms (DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (const edm::Run &r, const edm::EventSetup &c) override
 
 DrellYanValidation (const edm::ParameterSet &)
 
 ~DrellYanValidation () 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
 

Private Attributes

int _flavor
 decay flavor More...
 
std::string _name
 decay flavor name More...
 
MonitorElementcos_theta_gamma_lepton
 
MonitorElementdilep_mass
 
MonitorElementdilep_massPeak
 
MonitorElementdilep_pt
 
MonitorElementdilep_ptLog
 
MonitorElementdilep_rap
 
edm::ESHandle< HepPDT::ParticleDataTablefPDGTable
 PDT table. More...
 
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecordfPDGTableToken
 
MonitorElementgamma_energy
 
edm::InputTag hepmcCollection_
 
edm::EDGetTokenT< edm::HepMCProducthepmcCollectionToken_
 
MonitorElementleadeta
 
MonitorElementleadpt
 
MonitorElementnEvt
 
MonitorElementseceta
 
MonitorElementsecpt
 
WeightManager wmanager_
 
MonitorElementZdaughters
 
MonitorElementZmass
 
MonitorElementZmassPeak
 
MonitorElementZpt
 
MonitorElementZptLog
 
MonitorElementZrap
 

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 Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

Definition at line 32 of file DrellYanValidation.h.

Constructor & Destructor Documentation

◆ DrellYanValidation()

DrellYanValidation::DrellYanValidation ( const edm::ParameterSet iPSet)
explicit

Definition at line 19 of file DrellYanValidation.cc.

References fPDGTableToken, hepmcCollection_, and hepmcCollectionToken_.

20  : wmanager_(iPSet, consumesCollector()),
21  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
22  _flavor(iPSet.getParameter<int>("decaysTo")),
23  _name(iPSet.getParameter<std::string>("name")) {
24  hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
25  fPDGTableToken = esConsumes<edm::Transition::BeginRun>();
26 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > fPDGTableToken
edm::InputTag hepmcCollection_
int _flavor
decay flavor
WeightManager wmanager_
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
std::string _name
decay flavor name

◆ ~DrellYanValidation()

DrellYanValidation::~DrellYanValidation ( )
override

Definition at line 28 of file DrellYanValidation.cc.

28 {}

Member Function Documentation

◆ analyze()

void DrellYanValidation::analyze ( edm::Event const &  iEvent,
edm::EventSetup const &  iSetup 
)
overridevirtual

Gathering the HepMCProduct information

Reimplemented from DQMEDAnalyzer.

Definition at line 84 of file DrellYanValidation.cc.

References _flavor, funct::abs(), cms::cuda::assert(), funct::cos(), cos_theta_gamma_lepton, dilep_mass, dilep_massPeak, dilep_pt, dilep_ptLog, dilep_rap, dqm::impl::MonitorElement::Fill(), HepMCValidationHelper::findFSRPhotons(), fPDGTable, gamma_energy, edm::HepMCProduct::GetEvent(), hepmcCollectionToken_, mps_fire::i, iEvent, leadeta, leadpt, RazorAnalyzer::lep1, RazorAnalyzer::lep2, nEvt, LHEGenericFilter_cfi::ParticleID, HiggsValidation_cfi::pdg_id, perp(), edm::es::products(), seceta, secpt, jetUpdater_cfi::sort, HepMCValidationHelper::sortByPt(), submitPVValidationJobs::t, WeightManager::weight(), mps_merge::weight, wmanager_, x, y, z, Zdaughters, Zmass, ZmassPeak, Zpt, ZptLog, and Zrap.

84  {
85  // we *DO NOT* rely on a Z entry in the particle listings!
86 
89  iEvent.getByToken(hepmcCollectionToken_, evt);
90 
91  //Get EVENT
92  const HepMC::GenEvent* myGenEvent = evt->GetEvent();
93 
94  double weight = wmanager_.weight(iEvent);
95 
96  //std::cout << "weight: " << weight << std::endl;
97 
98  nEvt->Fill(0.5, weight);
99 
100  std::vector<const HepMC::GenParticle*> allproducts;
101 
102  //requires status 1 for leptons and neutrinos (except tau)
103  int requiredstatus = (std::abs(_flavor) == 11 || std::abs(_flavor) == 12 || std::abs(_flavor) == 13 ||
104  std::abs(_flavor) == 14 || std::abs(_flavor) == 16)
105  ? 1
106  : 3;
107 
108  bool vetotau =
109  true; //(std::abs(_flavor) == 11 || std::abs(_flavor) == 12 || std::abs(_flavor) ==13 || std::abs(_flavor) ==14 || std::abs(_flavor) ==16) ? true : false;
110 
111  for (HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin();
112  iter != myGenEvent->particles_end();
113  ++iter) {
114  if (vetotau) {
115  if ((*iter)->status() == 3 && std::abs((*iter)->pdg_id()) == 15)
116  return;
117  }
118  if ((*iter)->status() == requiredstatus) {
119  if (std::abs((*iter)->pdg_id()) == _flavor)
120  allproducts.push_back(*iter);
121  }
122  }
123 
124  //nothing to do if we don't have 2 particles
125  if (allproducts.size() < 2)
126  return;
127 
128  //sort them in pt
129  std::sort(allproducts.begin(), allproducts.end(), HepMCValidationHelper::sortByPt);
130 
131  //get the first element and the first following element with opposite charge
132  std::vector<const HepMC::GenParticle*> products;
133  products.push_back(allproducts.front());
134  const HepPDT::ParticleData* PData1 = fPDGTable->particle(HepPDT::ParticleID(allproducts.front()->pdg_id()));
135  double charge1 = PData1->charge();
136  for (unsigned int i = 1; i < allproducts.size(); ++i) {
137  const HepPDT::ParticleData* PData2 = fPDGTable->particle(HepPDT::ParticleID(allproducts[i]->pdg_id()));
138  double charge2 = PData2->charge();
139  if (charge1 * charge2 < 0)
140  products.push_back(allproducts[i]);
141  }
142 
143  //if we did not find any opposite charge pair there is nothing to do
144  if (products.size() < 2)
145  return;
146 
147  assert(products[0]->momentum().perp() >= products[1]->momentum().perp());
148 
149  //leading lepton with pt > 20.
150  if (products[0]->momentum().perp() < 20.)
151  return;
152 
153  //assemble FourMomenta
154  TLorentzVector lep1(products[0]->momentum().x(),
155  products[0]->momentum().y(),
156  products[0]->momentum().z(),
157  products[0]->momentum().t());
158  TLorentzVector lep2(products[1]->momentum().x(),
159  products[1]->momentum().y(),
160  products[1]->momentum().z(),
161  products[1]->momentum().t());
162  TLorentzVector dilepton_mom = lep1 + lep2;
163  TLorentzVector dilepton_andphoton_mom = dilepton_mom;
164 
165  //mass > 60.
166  if (dilepton_mom.M() < 60.)
167  return;
168 
169  //find possible qed fsr photons
170  std::vector<const HepMC::GenParticle*> fsrphotons;
171  HepMCValidationHelper::findFSRPhotons(products, myGenEvent, 0.1, fsrphotons);
172 
175 
176  std::vector<TLorentzVector> gammasMomenta;
177  for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho) {
178  TLorentzVector phomom(fsrphotons[ipho]->momentum().x(),
179  fsrphotons[ipho]->momentum().y(),
180  fsrphotons[ipho]->momentum().z(),
181  fsrphotons[ipho]->momentum().t());
182  dilepton_andphoton_mom += phomom;
183  Zdaughters->Fill(fsrphotons[ipho]->pdg_id(), weight);
184  gammasMomenta.push_back(phomom);
185  }
186  //Fill Z histograms
187  Zmass->Fill(dilepton_andphoton_mom.M(), weight);
188  ZmassPeak->Fill(dilepton_andphoton_mom.M(), weight);
189  Zpt->Fill(dilepton_andphoton_mom.Pt(), weight);
190  ZptLog->Fill(log10(dilepton_andphoton_mom.Pt()), weight);
191  Zrap->Fill(dilepton_andphoton_mom.Rapidity(), weight);
192 
193  //Fill dilepton histograms
194  dilep_mass->Fill(dilepton_mom.M(), weight);
195  dilep_massPeak->Fill(dilepton_mom.M(), weight);
196  dilep_pt->Fill(dilepton_mom.Pt(), weight);
197  dilep_ptLog->Fill(log10(dilepton_mom.Pt()), weight);
198  dilep_rap->Fill(dilepton_mom.Rapidity(), weight);
199 
200  //Fill lepton histograms
201  leadpt->Fill(lep1.Pt(), weight);
202  secpt->Fill(lep2.Pt(), weight);
203  leadeta->Fill(lep1.Eta(), weight);
204  seceta->Fill(lep2.Eta(), weight);
205 
206  //boost everything in the Z frame
207  TVector3 boost = dilepton_andphoton_mom.BoostVector();
208  boost *= -1.;
209  lep1.Boost(boost);
210  lep2.Boost(boost);
211  for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho) {
212  gammasMomenta[ipho].Boost(boost);
213  }
214  std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
215 
216  //fill gamma histograms
217  if (!gammasMomenta.empty() && dilepton_andphoton_mom.M() > 50.) {
218  gamma_energy->Fill(gammasMomenta.front().E(), weight);
219  double dphi = lep1.DeltaR(gammasMomenta.front()) < lep2.DeltaR(gammasMomenta.front())
220  ? lep1.DeltaPhi(gammasMomenta.front())
221  : lep2.DeltaPhi(gammasMomenta.front());
223  }
224 
225 } //analyze
MonitorElement * secpt
MonitorElement * gamma_energy
MonitorElement * dilep_rap
Definition: CLHEP.h:16
MonitorElement * Zmass
ESProducts< std::remove_reference_t< TArgs >... > products(TArgs &&... args)
Definition: ESProducts.h:128
MonitorElement * cos_theta_gamma_lepton
MonitorElement * nEvt
MonitorElement * leadeta
bool sortByPt(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
MonitorElement * Zrap
Definition: weight.py:1
assert(be >=bs)
MonitorElement * leadpt
void findFSRPhotons(const std::vector< const HepMC::GenParticle *> &leptons, const std::vector< const HepMC::GenParticle *> &all, double deltaR, std::vector< const HepMC::GenParticle *> &photons)
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
void Fill(long long x)
int _flavor
decay flavor
int iEvent
Definition: GenABIO.cc:224
MonitorElement * dilep_massPeak
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
WeightManager wmanager_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorElement * Zdaughters
MonitorElement * dilep_pt
HepPDT::ParticleData ParticleData
T perp() const
Magnitude of transverse component.
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
MonitorElement * ZptLog
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
MonitorElement * seceta
MonitorElement * dilep_mass
MonitorElement * ZmassPeak
MonitorElement * dilep_ptLog
lep1
print &#39;MRbb(1b)&#39;,event.mr_bb
double weight(const edm::Event &)
MonitorElement * Zpt

◆ bookHistograms()

void DrellYanValidation::bookHistograms ( DQMStore::IBooker i,
edm::Run const &  ,
edm::EventSetup const &   
)
overridevirtual

Setting the DQM top directories

Implements DQMEDAnalyzer.

Definition at line 34 of file DrellYanValidation.cc.

References _name, cos_theta_gamma_lepton, dilep_mass, dilep_massPeak, dilep_pt, dilep_ptLog, dilep_rap, ALCARECODTCalibSynchCosmicsDQM_cff::folderName, gamma_energy, mps_fire::i, leadeta, leadpt, nEvt, seceta, secpt, AlCaHLTBitMon_QueryRunRegistry::string, Zdaughters, Zmass, ZmassPeak, Zpt, ZptLog, and Zrap.

34  {
36  std::string folderName = "Generator/DrellYan";
37  folderName += _name;
38  DQMHelper dqm(&i);
39  i.setCurrentFolder(folderName);
40 
41  // Number of analyzed events
42  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1., "bin", "Number of Events");
43 
44  //Kinematics
45  Zmass = dqm.book1dHisto("Zmass", "inv. Mass Z", 70, 0, 140, "M_{Z} (GeV)", "Number of Events");
46  ZmassPeak = dqm.book1dHisto("ZmassPeak", "inv. Mass Z", 80, 80, 100, "M_{Z} (GeV)", "Number of Events");
47  Zpt = dqm.book1dHisto("Zpt", "Z pt", 100, 0, 200, "P_{t}^{Z} (GeV)", "Number of Events");
48  ZptLog =
49  dqm.book1dHisto("ZptLog", "log(Z pt)", 100, 0., 5., "log_{10}(P_{t}^{Z}) (log_{10}(GeV))", "Number of Events");
50  Zrap = dqm.book1dHisto("Zrap", "Z y", 100, -5, 5, "Y_{Z}", "Number of Events");
51  Zdaughters = dqm.book1dHisto("Zdaughters", "Z daughters", 60, -30, 30, "Z daughters (PDG ID)", "Number of Events");
52 
53  dilep_mass = dqm.book1dHisto("dilep_mass", "inv. Mass dilepton", 70, 0, 140, "M_{ll} (GeV)", "Number of Events");
55  dqm.book1dHisto("dilep_massPeak", "inv. Mass dilepton", 80, 80, 100, "M_{ll} (GeV)", "Number of Events");
56  dilep_pt = dqm.book1dHisto("dilep_pt", "dilepton pt", 100, 0, 200, "P_{t}^{ll} (GeV)", "Number of Events");
57  dilep_ptLog = dqm.book1dHisto(
58  "dilep_ptLog", "log(dilepton pt)", 100, 0., 5., "log_{10}(P_{t}^{ll}) (log_{10}(GeV))", "Number of Events");
59  dilep_rap = dqm.book1dHisto("dilep_rap", "dilepton y", 100, -5, 5, "Y_{ll}", "Number of Events");
60 
61  gamma_energy = dqm.book1dHisto("gamma_energy",
62  "photon energy in Z rest frame",
63  200,
64  0.,
65  100.,
66  "E_{#gamma}^{Z rest-frame} (GeV)",
67  "Number of Events");
68  cos_theta_gamma_lepton = dqm.book1dHisto("cos_theta_gamma_lepton",
69  "cos_theta_gamma_lepton in Z rest frame",
70  200,
71  -1,
72  1,
73  "cos(#theta_{#gamma-lepton}^{Z rest-frame})",
74  "Number of Events");
75 
76  leadpt = dqm.book1dHisto("leadpt", "leading lepton pt", 200, 0., 200., "P_{t}^{1st-lepton}", "Number of Events");
77  secpt = dqm.book1dHisto("secpt", "second lepton pt", 200, 0., 200., "P_{t}^{2nd-lepton}", "Number of Events");
78  leadeta = dqm.book1dHisto("leadeta", "leading lepton eta", 100, -5., 5., "#eta^{1st-lepton}", "Number of Events");
79  seceta = dqm.book1dHisto("seceta", "second lepton eta", 100, -5., 5., "#eta^{2nd-lepton}", "Number of Events");
80 
81  return;
82 }
MonitorElement * secpt
MonitorElement * gamma_energy
MonitorElement * dilep_rap
MonitorElement * Zmass
MonitorElement * cos_theta_gamma_lepton
MonitorElement * nEvt
MonitorElement * leadeta
MonitorElement * Zrap
MonitorElement * leadpt
MonitorElement * dilep_massPeak
MonitorElement * Zdaughters
MonitorElement * dilep_pt
MonitorElement * ZptLog
std::string _name
decay flavor name
MonitorElement * seceta
MonitorElement * dilep_mass
MonitorElement * ZmassPeak
MonitorElement * dilep_ptLog
MonitorElement * Zpt
Definition: DQMStore.h:18

◆ dqmBeginRun()

void DrellYanValidation::dqmBeginRun ( const edm::Run r,
const edm::EventSetup c 
)
overridevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 30 of file DrellYanValidation.cc.

References DummyCfis::c, fPDGTable, and fPDGTableToken.

30  {
31  fPDGTable = c.getHandle(fPDGTableToken);
32 }
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > fPDGTableToken
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.

Member Data Documentation

◆ _flavor

int DrellYanValidation::_flavor
private

decay flavor

Definition at line 56 of file DrellYanValidation.h.

Referenced by analyze().

◆ _name

std::string DrellYanValidation::_name
private

◆ cos_theta_gamma_lepton

MonitorElement * DrellYanValidation::cos_theta_gamma_lepton
private

Definition at line 53 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ dilep_mass

MonitorElement* DrellYanValidation::dilep_mass
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ dilep_massPeak

MonitorElement * DrellYanValidation::dilep_massPeak
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ dilep_pt

MonitorElement * DrellYanValidation::dilep_pt
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ dilep_ptLog

MonitorElement * DrellYanValidation::dilep_ptLog
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ dilep_rap

MonitorElement * DrellYanValidation::dilep_rap
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ fPDGTable

edm::ESHandle<HepPDT::ParticleDataTable> DrellYanValidation::fPDGTable
private

PDT table.

Definition at line 46 of file DrellYanValidation.h.

Referenced by analyze(), and dqmBeginRun().

◆ fPDGTableToken

edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> DrellYanValidation::fPDGTableToken
private

Definition at line 47 of file DrellYanValidation.h.

Referenced by dqmBeginRun(), and DrellYanValidation().

◆ gamma_energy

MonitorElement* DrellYanValidation::gamma_energy
private

Definition at line 53 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ hepmcCollection_

edm::InputTag DrellYanValidation::hepmcCollection_
private

Definition at line 43 of file DrellYanValidation.h.

Referenced by DrellYanValidation().

◆ hepmcCollectionToken_

edm::EDGetTokenT<edm::HepMCProduct> DrellYanValidation::hepmcCollectionToken_
private

Definition at line 60 of file DrellYanValidation.h.

Referenced by analyze(), and DrellYanValidation().

◆ leadeta

MonitorElement * DrellYanValidation::leadeta
private

Definition at line 52 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ leadpt

MonitorElement* DrellYanValidation::leadpt
private

Definition at line 52 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ nEvt

MonitorElement* DrellYanValidation::nEvt
private

Definition at line 49 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ seceta

MonitorElement * DrellYanValidation::seceta
private

Definition at line 52 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ secpt

MonitorElement * DrellYanValidation::secpt
private

Definition at line 52 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ wmanager_

WeightManager DrellYanValidation::wmanager_
private

Definition at line 42 of file DrellYanValidation.h.

Referenced by analyze().

◆ Zdaughters

MonitorElement * DrellYanValidation::Zdaughters
private

Definition at line 50 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ Zmass

MonitorElement* DrellYanValidation::Zmass
private

Definition at line 50 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ ZmassPeak

MonitorElement * DrellYanValidation::ZmassPeak
private

Definition at line 50 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ Zpt

MonitorElement * DrellYanValidation::Zpt
private

Definition at line 50 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ ZptLog

MonitorElement * DrellYanValidation::ZptLog
private

Definition at line 50 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

◆ Zrap

MonitorElement * DrellYanValidation::Zrap
private

Definition at line 50 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().