CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
DrellYanValidation Class Reference

#include <DrellYanValidation.h>

Inheritance diagram for DrellYanValidation:
DQMEDAnalyzer edm::stream::EDAnalyzer< edm::RunSummaryCache< dqmDetails::NoCache >, edm::LuminosityBlockSummaryCache< dqmDetails::NoCache > > edm::stream::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

virtual void analyze (edm::Event const &, edm::EventSetup const &) override
 
virtual void bookHistograms (DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
 
virtual void dqmBeginRun (const edm::Run &r, const edm::EventSetup &c) override
 
 DrellYanValidation (const edm::ParameterSet &)
 
virtual ~DrellYanValidation ()
 
- Public Member Functions inherited from DQMEDAnalyzer
virtual void beginRun (edm::Run const &, edm::EventSetup const &) final
 
virtual void beginStream (edm::StreamID id) final
 
 DQMEDAnalyzer (void)
 
virtual void endLuminosityBlockSummary (edm::LuminosityBlock const &, edm::EventSetup const &, dqmDetails::NoCache *) const final
 
virtual void endRunSummary (edm::Run const &, edm::EventSetup const &, dqmDetails::NoCache *) const final
 
uint32_t streamId () const
 
- Public Member Functions inherited from edm::stream::EDAnalyzer< edm::RunSummaryCache< dqmDetails::NoCache >, edm::LuminosityBlockSummaryCache< dqmDetails::NoCache > >
 EDAnalyzer ()=default
 
- Public Member Functions inherited from edm::stream::EDAnalyzerBase
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDAnalyzerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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::ParticleDataTable
fPDGTable
 PDT table. More...
 
MonitorElementgamma_energy
 
edm::InputTag hepmcCollection_
 
edm::EDGetTokenT
< edm::HepMCProduct
hepmcCollectionToken_
 
MonitorElementleadeta
 
MonitorElementleadpt
 
MonitorElementnEvt
 
MonitorElementseceta
 
MonitorElementsecpt
 
WeightManager wmanager_
 
MonitorElementZdaughters
 
MonitorElementZmass
 
MonitorElementZmassPeak
 
MonitorElementZpt
 
MonitorElementZptLog
 
MonitorElementZrap
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDAnalyzer< edm::RunSummaryCache< dqmDetails::NoCache >, edm::LuminosityBlockSummaryCache< dqmDetails::NoCache > >
typedef CacheContexts< T...> CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T...> HasAbility
 
typedef
CacheTypes::LuminosityBlockCache 
LuminosityBlockCache
 
typedef
LuminosityBlockContextT
< LuminosityBlockCache,
RunCache, GlobalCache
LuminosityBlockContext
 
typedef
CacheTypes::LuminosityBlockSummaryCache 
LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache,
GlobalCache
RunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Public Types inherited from edm::stream::EDAnalyzerBase
typedef EDAnalyzerAdaptorBase ModuleType
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static std::shared_ptr
< dqmDetails::NoCache
globalBeginLuminosityBlockSummary (edm::LuminosityBlock const &, edm::EventSetup const &, LuminosityBlockContext const *)
 
static std::shared_ptr
< dqmDetails::NoCache
globalBeginRunSummary (edm::Run const &, edm::EventSetup const &, RunContext const *)
 
static void globalEndLuminosityBlockSummary (edm::LuminosityBlock const &, edm::EventSetup const &, LuminosityBlockContext const *, dqmDetails::NoCache *)
 
static void globalEndRunSummary (edm::Run const &, edm::EventSetup const &, RunContext const *, dqmDetails::NoCache *)
 
- Static Public Member Functions inherited from edm::stream::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::stream::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 33 of file DrellYanValidation.h.

Constructor & Destructor Documentation

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

Definition at line 19 of file DrellYanValidation.cc.

References hepmcCollection_, and hepmcCollectionToken_.

19  :
21  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
22  _flavor(iPSet.getParameter<int>("decaysTo")),
23  _name(iPSet.getParameter<std::string>("name"))
24 {
25  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
26 }
T getParameter(std::string const &) const
edm::InputTag hepmcCollection_
int _flavor
decay flavor
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
WeightManager wmanager_
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
std::string _name
decay flavor name
DrellYanValidation::~DrellYanValidation ( )
virtual

Definition at line 28 of file DrellYanValidation.cc.

28 {}

Member Function Documentation

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

Gathering the HepMCProduct information

Implements edm::stream::EDAnalyzerBase.

Definition at line 69 of file DrellYanValidation.cc.

References _flavor, funct::abs(), funct::cos(), cos_theta_gamma_lepton, dilep_mass, dilep_massPeak, dilep_pt, dilep_ptLog, dilep_rap, MonitorElement::Fill(), HepMCValidationHelper::findFSRPhotons(), fPDGTable, gamma_energy, edm::Event::getByToken(), hepmcCollectionToken_, i, getDQMSummary::iter, leadeta, leadpt, nEvt, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, perp(), edm::es::products(), seceta, secpt, python.multivaluedict::sort(), HepMCValidationHelper::sortByPt(), edmStreamStallGrapher::t, WeightManager::weight(), histoStyle::weight, wmanager_, x, detailsBasic3DVector::y, detailsBasic3DVector::z, Zdaughters, Zmass, ZmassPeak, Zpt, ZptLog, and Zrap.

70 {
71 
72  // we *DO NOT* rely on a Z entry in the particle listings!
73 
76  iEvent.getByToken(hepmcCollectionToken_, evt);
77 
78  //Get EVENT
79  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
80 
81  double weight = wmanager_.weight(iEvent);
82 
83  //std::cout << "weight: " << weight << std::endl;
84 
85  nEvt->Fill(0.5,weight);
86 
87  std::vector<const HepMC::GenParticle*> allproducts;
88 
89  //requires status 1 for leptons and neutrinos (except tau)
90  int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? 1 : 3;
91 
92  bool vetotau = true; //(abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? true : false;
93 
94  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
95  if (vetotau) {
96  if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return;
97  }
98  if((*iter)->status()==requiredstatus) {
99  if(abs((*iter)->pdg_id())==_flavor)
100  allproducts.push_back(*iter);
101  }
102  }
103 
104  //nothing to do if we don't have 2 particles
105  if (allproducts.size() < 2) return;
106 
107  //sort them in pt
108  std::sort(allproducts.begin(), allproducts.end(), HepMCValidationHelper::sortByPt);
109 
110  //get the first element and the first following element with opposite charge
111  std::vector<const HepMC::GenParticle*> products;
112  products.push_back(allproducts.front());
113  const HepPDT::ParticleData* PData1 = fPDGTable->particle(HepPDT::ParticleID(allproducts.front()->pdg_id()));
114  double charge1 = PData1->charge();
115  for (unsigned int i = 1; i < allproducts.size(); ++i ){
116  const HepPDT::ParticleData* PData2 = fPDGTable->particle(HepPDT::ParticleID(allproducts[i]->pdg_id()));
117  double charge2 = PData2->charge();
118  if (charge1*charge2 < 0) products.push_back(allproducts[i]);
119  }
120 
121  //if we did not find any opposite charge pair there is nothing to do
122  if (products.size() < 2) return;
123 
124  assert(products[0]->momentum().perp() >= products[1]->momentum().perp());
125 
126  //leading lepton with pt > 20.
127  if (products[0]->momentum().perp() < 20.) return;
128 
129  //assemble FourMomenta
130  TLorentzVector lep1(products[0]->momentum().x(), products[0]->momentum().y(), products[0]->momentum().z(), products[0]->momentum().t());
131  TLorentzVector lep2(products[1]->momentum().x(), products[1]->momentum().y(), products[1]->momentum().z(), products[1]->momentum().t());
132  TLorentzVector dilepton_mom = lep1 + lep2;
133  TLorentzVector dilepton_andphoton_mom = dilepton_mom;
134 
135  //mass > 60.
136  if (dilepton_mom.M() < 60.) return;
137 
138  //find possible qed fsr photons
139  std::vector<const HepMC::GenParticle*> fsrphotons;
140  HepMCValidationHelper::findFSRPhotons(products, myGenEvent, 0.1, fsrphotons);
141 
142  Zdaughters->Fill(products[0]->pdg_id(),weight);
143  Zdaughters->Fill(products[1]->pdg_id(),weight);
144 
145  std::vector<TLorentzVector> gammasMomenta;
146  for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
147  TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t());
148  dilepton_andphoton_mom += phomom;
149  Zdaughters->Fill(fsrphotons[ipho]->pdg_id(),weight);
150  gammasMomenta.push_back(phomom);
151  }
152  //Fill Z histograms
153  Zmass->Fill(dilepton_andphoton_mom.M(),weight);
154  ZmassPeak->Fill(dilepton_andphoton_mom.M(),weight);
155  Zpt->Fill(dilepton_andphoton_mom.Pt(),weight);
156  ZptLog->Fill(log10(dilepton_andphoton_mom.Pt()),weight);
157  Zrap->Fill(dilepton_andphoton_mom.Rapidity(),weight);
158 
159  //Fill dilepton histograms
160  dilep_mass->Fill(dilepton_mom.M(),weight);
161  dilep_massPeak->Fill(dilepton_mom.M(),weight);
162  dilep_pt->Fill(dilepton_mom.Pt(),weight);
163  dilep_ptLog->Fill(log10(dilepton_mom.Pt()),weight);
164  dilep_rap->Fill(dilepton_mom.Rapidity(),weight);
165 
166  //Fill lepton histograms
167  leadpt->Fill(lep1.Pt(),weight);
168  secpt->Fill(lep2.Pt(),weight);
169  leadeta->Fill(lep1.Eta(),weight);
170  seceta->Fill(lep2.Eta(),weight);
171 
172  //boost everything in the Z frame
173  TVector3 boost = dilepton_andphoton_mom.BoostVector();
174  boost*=-1.;
175  lep1.Boost(boost);
176  lep2.Boost(boost);
177  for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
178  gammasMomenta[ipho].Boost(boost);
179  }
180  std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
181 
182  //fill gamma histograms
183  if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
184  gamma_energy->Fill(gammasMomenta.front().E(),weight);
185  double dphi = lep1.DeltaR(gammasMomenta.front()) < lep2.DeltaR(gammasMomenta.front()) ?
186  lep1.DeltaPhi(gammasMomenta.front()) : lep2.DeltaPhi(gammasMomenta.front());
187  cos_theta_gamma_lepton->Fill(cos(dphi),weight);
188  }
189 
190 }//analyze
MonitorElement * secpt
MonitorElement * gamma_energy
int i
Definition: DBlmapReader.cc:9
MonitorElement * dilep_rap
void findFSRPhotons(const std::vector< const HepMC::GenParticle * > &leptons, const std::vector< const HepMC::GenParticle * > &all, double deltaR, std::vector< const HepMC::GenParticle * > &photons)
MonitorElement * Zmass
MonitorElement * cos_theta_gamma_lepton
MonitorElement * nEvt
MonitorElement * leadeta
bool sortByPt(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
MonitorElement * Zrap
MonitorElement * leadpt
float float float z
void Fill(long long x)
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
int _flavor
decay flavor
int iEvent
Definition: GenABIO.cc:230
MonitorElement * dilep_massPeak
ESProducts< T, S > products(const T &i1, const S &i2)
Definition: ESProducts.h:189
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
MonitorElement * ZptLog
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
MonitorElement * seceta
T perp() const
Magnitude of transverse component.
MonitorElement * dilep_mass
MonitorElement * ZmassPeak
MonitorElement * dilep_ptLog
Definition: DDAxes.h:10
int weight
Definition: histoStyle.py:50
double weight(const edm::Event &)
MonitorElement * Zpt
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, DQMStore::IBooker::book1D(), cos_theta_gamma_lepton, dilep_mass, dilep_massPeak, dilep_pt, dilep_ptLog, dilep_rap, gamma_energy, leadeta, leadpt, nEvt, seceta, secpt, DQMStore::IBooker::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, Zdaughters, Zmass, ZmassPeak, Zpt, ZptLog, and Zrap.

34  {
35 
37  std::string folderName = "Generator/DrellYan";
38  folderName+=_name;
39  i.setCurrentFolder(folderName.c_str());
40 
41  // Number of analyzed events
42  nEvt = i.book1D("nEvt", "n analyzed Events", 1, 0., 1.);
43 
44  //Kinematics
45  Zmass = i.book1D("Zmass","inv. Mass Z", 70 ,0,140);
46  ZmassPeak = i.book1D("ZmassPeak","inv. Mass Z", 80 ,80 ,100);
47  Zpt = i.book1D("Zpt","Z pt",100,0,200);
48  ZptLog = i.book1D("ZptLog","log(Z pt)",100,0.,5.);
49  Zrap = i.book1D("Zrap", "Z y", 100, -5, 5);
50  Zdaughters = i.book1D("Zdaughters", "Z daughters", 60, -30, 30);
51 
52  dilep_mass = i.book1D("dilep_mass","inv. Mass dilepton", 70 ,0,140);
53  dilep_massPeak = i.book1D("dilep_massPeak","inv. Mass dilepton", 80 ,80 ,100);
54  dilep_pt = i.book1D("dilep_pt","dilepton pt",100,0,200);
55  dilep_ptLog = i.book1D("dilep_ptLog","log(dilepton pt)",100,0.,5.);
56  dilep_rap = i.book1D("dilep_rap", "dilepton y", 100, -5, 5);
57 
58  gamma_energy = i.book1D("gamma_energy", "photon energy in Z rest frame", 200, 0., 100.);
59  cos_theta_gamma_lepton = i.book1D("cos_theta_gamma_lepton", "cos_theta_gamma_lepton in Z rest frame", 200, -1, 1);
60 
61  leadpt = i.book1D("leadpt","leading lepton pt", 200, 0., 200.);
62  secpt = i.book1D("secpt","second lepton pt", 200, 0., 200.);
63  leadeta = i.book1D("leadeta","leading lepton eta", 100, -5., 5.);
64  seceta = i.book1D("seceta","second lepton eta", 100, -5., 5.);
65 
66  return;
67 }
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 * book1D(Args &&...args)
Definition: DQMStore.h:113
MonitorElement * Zdaughters
MonitorElement * dilep_pt
MonitorElement * ZptLog
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
std::string _name
decay flavor name
MonitorElement * seceta
MonitorElement * dilep_mass
MonitorElement * ZmassPeak
MonitorElement * dilep_ptLog
MonitorElement * Zpt
void DrellYanValidation::dqmBeginRun ( const edm::Run r,
const edm::EventSetup c 
)
overridevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 30 of file DrellYanValidation.cc.

References fPDGTable, and edm::EventSetup::getData().

30  {
31  c.getData( fPDGTable );
32 }
void getData(T &iHolder) const
Definition: EventSetup.h:67
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.

Member Data Documentation

int DrellYanValidation::_flavor
private

decay flavor

Definition at line 57 of file DrellYanValidation.h.

Referenced by analyze().

std::string DrellYanValidation::_name
private

decay flavor name

Definition at line 59 of file DrellYanValidation.h.

Referenced by bookHistograms().

MonitorElement * DrellYanValidation::cos_theta_gamma_lepton
private

Definition at line 54 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* DrellYanValidation::dilep_mass
private

Definition at line 52 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::dilep_massPeak
private

Definition at line 52 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::dilep_pt
private

Definition at line 52 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::dilep_ptLog
private

Definition at line 52 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::dilep_rap
private

Definition at line 52 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

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

PDT table.

Definition at line 48 of file DrellYanValidation.h.

Referenced by analyze(), and dqmBeginRun().

MonitorElement* DrellYanValidation::gamma_energy
private

Definition at line 54 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

edm::InputTag DrellYanValidation::hepmcCollection_
private

Definition at line 45 of file DrellYanValidation.h.

Referenced by DrellYanValidation().

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

Definition at line 61 of file DrellYanValidation.h.

Referenced by analyze(), and DrellYanValidation().

MonitorElement * DrellYanValidation::leadeta
private

Definition at line 53 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* DrellYanValidation::leadpt
private

Definition at line 53 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* DrellYanValidation::nEvt
private

Definition at line 50 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::seceta
private

Definition at line 53 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::secpt
private

Definition at line 53 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

WeightManager DrellYanValidation::wmanager_
private

Definition at line 44 of file DrellYanValidation.h.

Referenced by analyze().

MonitorElement * DrellYanValidation::Zdaughters
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* DrellYanValidation::Zmass
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::ZmassPeak
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::Zpt
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::ZptLog
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * DrellYanValidation::Zrap
private

Definition at line 51 of file DrellYanValidation.h.

Referenced by analyze(), and bookHistograms().