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 Member Functions | Private Attributes
ZMCHistogrammer Class Reference
Inheritance diagram for ZMCHistogrammer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 ZMCHistogrammer (const edm::ParameterSet &pset)
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 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
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

virtual void analyze (const edm::Event &event, const edm::EventSetup &setup) override
 

Private Attributes

double angMax_
 
edm::EDGetTokenT
< reco::CandidateView
genToken_
 
TH1F * h_etaResZ_
 
TH1F * h_etaZ_
 
TH1F * h_etaZMC_
 
TH1F * h_invmMuMu_
 
TH1F * h_invmMuMuMC_
 
TH1F * h_mRatioZMuMu_
 
TH1F * h_mRatioZMuMuMC_
 
TH1F * h_mResZ_
 
TH1F * h_mResZMuMu_
 
TH1F * h_mResZMuMuMC_
 
TH1F * h_mZ_
 
TH1F * h_mZMC_
 
TH1F * h_nZ_
 
TH1F * h_nZMC_
 
TH1F * h_phiResZ_
 
TH1F * h_phiZ_
 
TH1F * h_phiZMC_
 
TH1F * h_ptResZ_
 
TH1F * h_ptZ_
 
TH1F * h_ptZMC_
 
TH1F * h_rapidityResZ_
 
TH1F * h_rapidityZ_
 
TH1F * h_rapidityZMC_
 
TH1F * h_thetaResZ_
 
TH1F * h_thetaZ_
 
TH1F * h_thetaZMC_
 
TH1F * h_weight_histo
 
edm::EDGetTokenT
< edm::HepMCProduct
hepMCProductToken_
 
bool isMCatNLO_
 
double massMax_
 
double massResMax_
 
edm::EDGetTokenT< std::vector
< reco::GenParticleRef > > 
matchToken_
 
unsigned int nbinsAng_
 
unsigned int nbinsMass_
 
unsigned int nbinsMassRes_
 
unsigned int nbinsPt_
 
double ptMax_
 
edm::EDGetTokenT
< reco::CandidateView
zToken_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- 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 10 of file ZMCHistogrammer.cc.

Constructor & Destructor Documentation

ZMCHistogrammer::ZMCHistogrammer ( const edm::ParameterSet pset)

Definition at line 48 of file ZMCHistogrammer.cc.

References angMax_, gather_cfg::cout, h_etaResZ_, h_etaZ_, h_etaZMC_, h_invmMuMu_, h_invmMuMuMC_, h_mRatioZMuMu_, h_mRatioZMuMuMC_, h_mResZ_, h_mResZMuMu_, h_mResZMuMuMC_, h_mZ_, h_mZMC_, h_nZ_, h_nZMC_, h_phiResZ_, h_phiZ_, h_phiZMC_, h_ptResZ_, h_ptZ_, h_ptZMC_, h_rapidityResZ_, h_rapidityZ_, h_rapidityZMC_, h_thetaResZ_, h_thetaZ_, h_thetaZMC_, h_weight_histo, TFileDirectory::make(), massMax_, massResMax_, TFileService::mkdir(), nbinsAng_, nbinsMass_, nbinsMassRes_, nbinsPt_, and ptMax_.

48  :
49  zToken_(consumes<CandidateView>(pset.getParameter<InputTag>("z"))),
50  genToken_(consumes<CandidateView>(pset.getParameter<InputTag>("gen"))),
51  matchToken_(consumes<std::vector<GenParticleRef> >(pset.getParameter<InputTag>("match"))),
52  hepMCProductToken_(consumes<HepMCProduct>(pset.getParameter<InputTag>("hepMCProductTag"))),
53  nbinsMass_(pset.getUntrackedParameter<unsigned int>("nbinsMass")),
54  nbinsPt_(pset.getUntrackedParameter<unsigned int>("nbinsPt")),
55  nbinsAng_(pset.getUntrackedParameter<unsigned int>("nbinsAng")),
56  nbinsMassRes_(pset.getUntrackedParameter<unsigned int>("nbinsMassRes")),
57  massMax_(pset.getUntrackedParameter<double>("massMax")),
58  ptMax_(pset.getUntrackedParameter<double>("ptMax")),
59  angMax_(pset.getUntrackedParameter<double>("angMax")),
60  massResMax_(pset.getUntrackedParameter<double>("massResMax")),
61  isMCatNLO_(pset.getUntrackedParameter<bool>("isMCatNLO")) {
62  cout << ">>> Z Histogrammer constructor" << endl;
64  TFileDirectory ZHisto = fs->mkdir( "ZRecoHisto" );
65  TFileDirectory ZMCHisto = fs->mkdir( "ZMCHisto" );
66  TFileDirectory ZResHisto = fs->mkdir( "ZResolutionHisto" );
67  //TFileDirectory Z2vs3MCHisto = fs->mkdir( "Z2vs3MCHisto" );
68  h_nZ_ = ZHisto.make<TH1F>("ZNumber", "number of Z particles", 11, -0.5, 10.5);
69  h_mZ_ = ZHisto.make<TH1F>("ZMass", "Z mass (GeV/c^{2})", nbinsMass_, 0, massMax_);
70  h_ptZ_ = ZHisto.make<TH1F>("ZPt", "Z p_{t} (GeV/c)", nbinsPt_, 0, ptMax_);
71  h_phiZ_ = ZHisto.make<TH1F>("ZPhi", "Z #phi", nbinsAng_, -angMax_, angMax_);
72  h_thetaZ_ = ZHisto.make<TH1F>("Ztheta", "Z #theta", nbinsAng_, 0, angMax_);
73  h_etaZ_ = ZHisto.make<TH1F>("ZEta", "Z #eta", nbinsAng_, -angMax_, angMax_);
74  h_rapidityZ_ = ZHisto.make<TH1F>("ZRapidity", "Z rapidity", nbinsAng_, -angMax_, angMax_);
75  h_invmMuMu_ = ZHisto.make<TH1F>("MuMuMass", "#mu #mu invariant mass",nbinsMass_, 0, massMax_);
76  h_weight_histo = ZHisto.make<TH1F>("weight_histo","weight_histo",20,-10,10);
77  h_nZMC_ = ZMCHisto.make<TH1F>("ZMCNumber", "number of Z MC particles", 11, -0.5, 10.5);
78  h_mZMC_ = ZMCHisto.make<TH1F>("ZMCMass", "Z MC mass (GeV/c^{2})", nbinsMass_, 0, massMax_);
79  h_ptZMC_ = ZMCHisto.make<TH1F>("ZMCPt", "Z MC p_{t} (GeV/c)", nbinsPt_, 0, ptMax_);
80  h_phiZMC_ = ZMCHisto.make<TH1F>("ZMCPhi", "Z MC #phi", nbinsAng_, -angMax_, angMax_);
81  h_thetaZMC_ = ZMCHisto.make<TH1F>("ZMCTheta", "Z MC #theta", nbinsAng_, 0, angMax_);
82  h_etaZMC_ = ZMCHisto.make<TH1F>("ZMCEta", "Z MC #eta", nbinsAng_, -angMax_, angMax_);
83  h_rapidityZMC_ = ZMCHisto.make<TH1F>("ZMCRapidity", "Z MC rapidity", nbinsAng_, -angMax_, angMax_);
84  h_invmMuMuMC_ = ZMCHisto.make<TH1F>("MuMuMCMass", "#mu #mu MC invariant mass",
85  nbinsMass_, 0, massMax_);
86  /*
87  h_mZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCMass", "Z MC st 2 vs st 3 mass (GeV/c^{2})",
88  nbinsMassRes_, -massResMax_, massResMax_);
89  h_ptZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCPt", "Z MC st 2 vs st 3 p_{t} (GeV/c)",
90  nbinsPt_, -ptMax_, ptMax_);
91  h_phiZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCPhi", "Z MC st 2 vs st 3 #phi",
92  nbinsAng_, -angMax_, angMax_);
93  h_thetaZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCTheta", "Z MC st 2 vs st 3 #theta",
94  nbinsAng_, -angMax_, angMax_);
95  h_etaZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCEta", "Z MC st 2 vs st 3 #eta",
96  nbinsAng_, -angMax_, angMax_);
97  h_rapidityZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCRapidity", "Z MC st 2 vs st 3 rapidity",
98  nbinsAng_, -angMax_, angMax_);
99  */
100  h_mResZ_ = ZResHisto.make<TH1F>("ZMassResolution", "Z mass Resolution (GeV/c^{2})",
102  h_ptResZ_ = ZResHisto.make<TH1F>("ZPtResolution", "Z p_{t} Resolution (GeV/c)",
103  nbinsPt_, -ptMax_, ptMax_);
104  h_phiResZ_ = ZResHisto.make<TH1F>("ZPhiResolution", "Z #phi Resolution",
106  h_thetaResZ_ = ZResHisto.make<TH1F>("ZThetaResolution", "Z #theta Resolution",
108  h_etaResZ_ = ZResHisto.make<TH1F>("ZEtaResolution", "Z #eta Resolution",
110  h_rapidityResZ_ = ZResHisto.make<TH1F>("ZRapidityResolution", "Z rapidity Resolution",
112  h_mResZMuMu_ = ZResHisto.make<TH1F>("ZToMuMuRecoMassResolution",
113  "Z Reco vs matched final state #mu #mu mass Difference (GeV/c^{2})",
115  h_mRatioZMuMu_ = ZResHisto.make<TH1F>("ZToMuMuRecoMassRatio",
116  "Z Reco vs matched final state #mu #mu mass Ratio",
117  4000, 0, 2);
118  h_mResZMuMuMC_ = ZResHisto.make<TH1F>("ZToMuMuMCMassResolution",
119  "Z vs final state #mu #mu MC mass Difference (GeV/c^{2})",
120  nbinsMassRes_/2 + 1, -2*massResMax_/nbinsMassRes_, massResMax_);
121  h_mRatioZMuMuMC_ = ZResHisto.make<TH1F>("ZToMuMuMCMassRatio",
122  "Z vs final state #mu #mu MC mass Ratio",
123  2002, 0.999, 2);
124 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< std::vector< reco::GenParticleRef > > matchToken_
edm::EDGetTokenT< reco::CandidateView > zToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
unsigned int nbinsAng_
T * make(const Args &...args) const
make new ROOT object
unsigned int nbinsMass_
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
edm::EDGetTokenT< reco::CandidateView > genToken_
tuple cout
Definition: gather_cfg.py:121
unsigned int nbinsPt_
edm::EDGetTokenT< edm::HepMCProduct > hepMCProductToken_
unsigned int nbinsMassRes_

Member Function Documentation

void ZMCHistogrammer::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
overrideprivatevirtual

Implements edm::EDAnalyzer.

Definition at line 126 of file ZMCHistogrammer.cc.

References funct::abs(), assert(), reco::Candidate::charge(), gather_cfg::cout, reco::Candidate::daughter(), reco::Candidate::eta(), relval_steps::gen(), genToken_, h_etaResZ_, h_etaZ_, h_etaZMC_, h_invmMuMu_, h_invmMuMuMC_, h_mRatioZMuMu_, h_mRatioZMuMuMC_, h_mResZ_, h_mResZMuMu_, h_mResZMuMuMC_, h_mZ_, h_mZMC_, h_nZ_, h_nZMC_, h_phiResZ_, h_phiZ_, h_phiZMC_, h_ptResZ_, h_ptZ_, h_ptZMC_, h_rapidityResZ_, h_rapidityZ_, h_rapidityZMC_, h_thetaResZ_, h_thetaZ_, h_thetaZMC_, h_weight_histo, hepMCProductToken_, i, isMCatNLO_, edm::Ref< C, T, F >::isNonnull(), j, relval_steps::k, reco::Candidate::mass(), match(), matchToken_, reco::Candidate::numberOfDaughters(), reco::Candidate::p4(), reco::Candidate::pdgId(), reco::Candidate::phi(), reco::Candidate::pt(), reco::Candidate::rapidity(), reco::Candidate::status(), reco::Candidate::theta(), histoStyle::weight, create_public_pileup_plots::weights, detailsBasic3DVector::z, and zToken_.

126  {
127  cout << ">>> Z Histogrammer analyze" << endl;
131  event.getByToken(zToken_, z);
132  event.getByToken(genToken_, gen);
133  event.getByToken(matchToken_, match);
134  h_nZ_->Fill(z->size());
135 
136 // get HepMC::GenEvent ...
137  Handle<HepMCProduct> evt_h;
138  event.getByToken(hepMCProductToken_, evt_h);
139  HepMC::GenEvent* evt = new HepMC::GenEvent(*(evt_h->GetEvent()));
140 
141 
142  // get weight and fill it to histogram
143  HepMC::WeightContainer weights = evt->weights();
144  double weight = weights.front();
145  if(!weight) weight=1.;
146  h_weight_histo->Fill(weight);
147 
148  if(isMCatNLO_) {
149  weight > 0 ? weight=1. : weight=-1.;
150  }
151 
152  for(unsigned int i = 0; i < z->size(); ++i) {
153  const Candidate &zCand = (*z)[i];
154  h_mZ_->Fill(zCand.mass(),weight );
155  h_ptZ_->Fill(zCand.pt(),weight);
156  h_phiZ_->Fill(zCand.phi(),weight);
157  h_thetaZ_->Fill(zCand.theta(),weight);
158  h_etaZ_->Fill(zCand.eta(),weight);
159  h_rapidityZ_->Fill(zCand.rapidity(),weight);
160  CandidateBaseRef zCandRef = z->refAt(i);
161 
162  GenParticleRef zMCMatch = (*match)[i];
163  if(zMCMatch.isNonnull() && zMCMatch->pdgId()==23) {
164  h_mResZ_->Fill(zCandRef->mass() - zMCMatch->mass());
165  h_ptResZ_->Fill(zCandRef->pt() - zMCMatch->pt());
166  h_phiResZ_->Fill(zCandRef->phi() - zMCMatch->phi());
167  h_thetaResZ_->Fill(zCandRef->theta() - zMCMatch->theta());
168  h_etaResZ_->Fill(zCandRef->eta() - zMCMatch->eta());
169  h_rapidityResZ_->Fill(zCandRef->rapidity() - zMCMatch->rapidity());
170  const Candidate * dau0 = zMCMatch->daughter(0);
171  const Candidate * dau1 = zMCMatch->daughter(1);
172  for(unsigned int i0 = 0; i0 < dau0->numberOfDaughters(); ++i0) {
173  const Candidate * ddau0 = dau0->daughter(i0);
174  if(abs(ddau0->pdgId())==13 && ddau0->status()==1) {
175  dau0 = ddau0; break;
176  }
177  }
178  for(unsigned int i1 = 0; i1 < dau1->numberOfDaughters(); ++i1) {
179  const Candidate * ddau1 = dau1->daughter(i1);
180  if(abs(ddau1->pdgId())==13 && ddau1->status()==1) {
181  dau1 = ddau1; break;
182  }
183  }
184  assert(abs(dau0->pdgId())==13 && dau0->status()==1);
185  assert(abs(dau1->pdgId())==13 && dau1->status()==1);
186  double invMass = (dau0->p4()+dau1->p4()).mass();
187  h_invmMuMu_->Fill(invMass,weight);
188  h_mResZMuMu_->Fill(zCand.mass() - invMass);
189  h_mRatioZMuMu_->Fill(zCand.mass()/invMass);
190  }
191  }
192  h_nZMC_->Fill(gen->size());
193  for(unsigned int i = 0; i < gen->size(); ++i) {
194  const Candidate &genCand = (*gen)[i];
195  if((genCand.pdgId() == 23) && (genCand.status() == 2)) //this is an intermediate Z0
196  cout << ">>> intermediate Z0 found, with " << genCand.numberOfDaughters()
197  << " daughters" << endl;
198  if((genCand.pdgId() == 23)&&(genCand.status() == 3)) { //this is a Z0
199  cout << ">>> Z0 found, with " << genCand.numberOfDaughters()
200  << " daughters" << endl;
201  h_mZMC_->Fill(genCand.mass(),weight);
202  h_ptZMC_->Fill(genCand.pt(),weight);
203  h_phiZMC_->Fill(genCand.phi(),weight);
204  h_thetaZMC_->Fill(genCand.theta(),weight);
205  h_etaZMC_->Fill(genCand.eta(),weight);
206  h_rapidityZMC_->Fill(genCand.rapidity(),weight);
207  Particle::LorentzVector pZ(0, 0, 0, 0);
208  int nMu = 0;
209  for(unsigned int j = 0; j < genCand.numberOfDaughters(); ++j) {
210  const Candidate *dauGen = genCand.daughter(j);
211  /*
212  if((dauGen->pdgId() == 23) && (dauGen->status() == 2)) {
213  h_mZ2vs3MC_->Fill(genCand.mass() - dauGen->mass());
214  h_ptZ2vs3MC_->Fill(genCand.pt() - dauGen->pt());
215  h_phiZ2vs3MC_->Fill(genCand.phi() - dauGen->phi());
216  h_thetaZ2vs3MC_->Fill(genCand.theta() - dauGen->theta());
217  h_etaZ2vs3MC_->Fill(genCand.eta() - dauGen->eta());
218  h_rapidityZ2vs3MC_->Fill(genCand.rapidity() - dauGen->rapidity());
219  }
220  */
221  if((abs(dauGen->pdgId()) == 13) && (dauGen->numberOfDaughters() != 0)) {
222  //we are looking for photons of final state radiation
223  cout << ">>> The muon " << j
224  << " has " << dauGen->numberOfDaughters() << " daughters" <<endl;
225  for(unsigned int k = 0; k < dauGen->numberOfDaughters(); ++k) {
226  const Candidate * dauMuGen = dauGen->daughter(k);
227  cout << ">>> Mu " << j
228  << " daughter MC " << k
229  << " PDG Id " << dauMuGen->pdgId()
230  << ", status " << dauMuGen->status()
231  << ", charge " << dauMuGen->charge()
232  << endl;
233  if(abs(dauMuGen->pdgId()) == 13 && dauMuGen->status() ==1) {
234  pZ += dauMuGen->p4();
235  nMu ++;
236  }
237  }
238  }
239  }
240  assert(nMu == 2);
241  double mZ = pZ.mass();
242  h_invmMuMuMC_->Fill(mZ,weight);
243  h_mResZMuMuMC_->Fill(genCand.mass() - mZ);
244  h_mRatioZMuMuMC_->Fill(genCand.mass()/mZ);
245  }
246  }
247 }
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
virtual double pt() const =0
transverse momentum
virtual double mass() const =0
mass
virtual int status() const =0
status word
assert(m_qm.get())
virtual double rapidity() const =0
rapidity
edm::EDGetTokenT< std::vector< reco::GenParticleRef > > matchToken_
edm::EDGetTokenT< reco::CandidateView > zToken_
float float float z
virtual size_type numberOfDaughters() const =0
number of daughters
virtual double theta() const =0
momentum polar angle
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
virtual int charge() const =0
electric charge
virtual int pdgId() const =0
PDG identifier.
edm::EDGetTokenT< reco::CandidateView > genToken_
tuple cout
Definition: gather_cfg.py:121
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
int weight
Definition: histoStyle.py:50
edm::EDGetTokenT< edm::HepMCProduct > hepMCProductToken_
math::PtEtaPhiELorentzVectorF LorentzVector
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector

Member Data Documentation

double ZMCHistogrammer::angMax_
private

Definition at line 20 of file ZMCHistogrammer.cc.

Referenced by ZMCHistogrammer().

edm::EDGetTokenT<reco::CandidateView> ZMCHistogrammer::genToken_
private

Definition at line 16 of file ZMCHistogrammer.cc.

Referenced by analyze().

TH1F * ZMCHistogrammer::h_etaResZ_
private

Definition at line 26 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_etaZ_
private

Definition at line 21 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_etaZMC_
private

Definition at line 23 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F* ZMCHistogrammer::h_invmMuMu_
private

Definition at line 22 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F* ZMCHistogrammer::h_invmMuMuMC_
private

Definition at line 24 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_mRatioZMuMu_
private

Definition at line 27 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_mRatioZMuMuMC_
private

Definition at line 28 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F* ZMCHistogrammer::h_mResZ_
private

Definition at line 26 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F* ZMCHistogrammer::h_mResZMuMu_
private

Definition at line 27 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F* ZMCHistogrammer::h_mResZMuMuMC_
private

Definition at line 28 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_mZ_
private

Definition at line 21 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_mZMC_
private

Definition at line 23 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F* ZMCHistogrammer::h_nZ_
private

Definition at line 21 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F* ZMCHistogrammer::h_nZMC_
private

Definition at line 23 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_phiResZ_
private

Definition at line 26 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_phiZ_
private

Definition at line 21 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_phiZMC_
private

Definition at line 23 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_ptResZ_
private

Definition at line 26 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_ptZ_
private

Definition at line 21 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_ptZMC_
private

Definition at line 23 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_rapidityResZ_
private

Definition at line 26 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_rapidityZ_
private

Definition at line 21 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_rapidityZMC_
private

Definition at line 23 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_thetaResZ_
private

Definition at line 26 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_thetaZ_
private

Definition at line 21 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_thetaZMC_
private

Definition at line 23 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

TH1F * ZMCHistogrammer::h_weight_histo
private

Definition at line 24 of file ZMCHistogrammer.cc.

Referenced by analyze(), and ZMCHistogrammer().

edm::EDGetTokenT<edm::HepMCProduct> ZMCHistogrammer::hepMCProductToken_
private

Definition at line 18 of file ZMCHistogrammer.cc.

Referenced by analyze().

bool ZMCHistogrammer::isMCatNLO_
private

Definition at line 29 of file ZMCHistogrammer.cc.

Referenced by analyze().

double ZMCHistogrammer::massMax_
private

Definition at line 20 of file ZMCHistogrammer.cc.

Referenced by ZMCHistogrammer().

double ZMCHistogrammer::massResMax_
private

Definition at line 20 of file ZMCHistogrammer.cc.

Referenced by ZMCHistogrammer().

edm::EDGetTokenT<std::vector<reco::GenParticleRef> > ZMCHistogrammer::matchToken_
private

Definition at line 17 of file ZMCHistogrammer.cc.

Referenced by analyze().

unsigned int ZMCHistogrammer::nbinsAng_
private

Definition at line 19 of file ZMCHistogrammer.cc.

Referenced by ZMCHistogrammer().

unsigned int ZMCHistogrammer::nbinsMass_
private

Definition at line 19 of file ZMCHistogrammer.cc.

Referenced by ZMCHistogrammer().

unsigned int ZMCHistogrammer::nbinsMassRes_
private

Definition at line 19 of file ZMCHistogrammer.cc.

Referenced by ZMCHistogrammer().

unsigned int ZMCHistogrammer::nbinsPt_
private

Definition at line 19 of file ZMCHistogrammer.cc.

Referenced by ZMCHistogrammer().

double ZMCHistogrammer::ptMax_
private

Definition at line 20 of file ZMCHistogrammer.cc.

Referenced by ZMCHistogrammer().

edm::EDGetTokenT<reco::CandidateView> ZMCHistogrammer::zToken_
private

Definition at line 15 of file ZMCHistogrammer.cc.

Referenced by analyze().