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
gamma_radiative_analyzer Class Reference
Inheritance diagram for gamma_radiative_analyzer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 gamma_radiative_analyzer (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 ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector
< ProductResolverIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

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

Private Attributes

double dRTrk_
 
double dRVeto_
 
bool FSR_mu
 
bool FSR_mu0
 
bool FSR_mu1
 
bool FSR_tk
 
TH2D * h_gamma_pt_eta_
 
TH2D * h_mu_pt_eta_FSR_
 
TH2D * h_mu_pt_eta_no_FSR_
 
int numOfEvent
 
int numofGamma
 
double ptThreshold_
 
int zmmcounter
 
int zmscounter
 
int zmtcounter
 
EDGetTokenT< GenParticleMatchzMuMuMatchMapToken_
 
EDGetTokenT< CandidateViewzMuMuToken_
 
EDGetTokenT< GenParticleMatchzMuSaMatchMapToken_
 
EDGetTokenT< CandidateViewzMuSaToken_
 
EDGetTokenT< GenParticleMatchzMuTkMatchMapToken_
 
EDGetTokenT< CandidateViewzMuTkToken_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- 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 52 of file gamma_radiative_analysis.cc.

Constructor & Destructor Documentation

gamma_radiative_analyzer::gamma_radiative_analyzer ( const edm::ParameterSet pset)

Definition at line 76 of file gamma_radiative_analysis.cc.

References fs, h_gamma_pt_eta_, h_mu_pt_eta_FSR_, h_mu_pt_eta_no_FSR_, TFileService::make(), numOfEvent, numofGamma, zmmcounter, zmscounter, and zmtcounter.

76  :
77  zMuMuToken_(consumes<CandidateView>(pset.getParameter<InputTag>("zMuMu"))),
78  zMuMuMatchMapToken_(mayConsume<GenParticleMatch>(pset.getParameter<InputTag>("zMuMuMatchMap"))),
79  zMuTkToken_(consumes<CandidateView>(pset.getParameter<InputTag>("zMuTk"))),
80  zMuTkMatchMapToken_(mayConsume<GenParticleMatch>(pset.getParameter<InputTag>("zMuTkMatchMap"))),
81  zMuSaToken_(consumes<CandidateView>(pset.getParameter<InputTag>("zMuSa"))),
82  zMuSaMatchMapToken_(mayConsume<GenParticleMatch>(pset.getParameter<InputTag>("zMuSaMatchMap"))){
83  zmmcounter=0;
84  zmscounter=0;
85  zmtcounter=0;
86  numOfEvent=0;
87  numofGamma=0;
89 
90  // general histograms
91  h_gamma_pt_eta_= fs->make<TH2D>("h_gamma_pt_eta","pt vs eta of gamma",100,20,100,100,-2.0,2.0);
92  h_mu_pt_eta_FSR_= fs->make<TH2D>("h_mu_pt_eta_FSR","pt vs eta of muon with FSR",100,20,100,100,-2.0,2.0 );
93  h_mu_pt_eta_no_FSR_= fs->make<TH2D>("h_mu_pt_eta_no_FSR","pt vs eta of of muon withot FSR",100,20,100,100,-2.0,2.0);
94 }
T getParameter(std::string const &) const
EDGetTokenT< CandidateView > zMuMuToken_
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
EDGetTokenT< CandidateView > zMuSaToken_
EDGetTokenT< GenParticleMatch > zMuMuMatchMapToken_
EDGetTokenT< GenParticleMatch > zMuTkMatchMapToken_
EDGetTokenT< GenParticleMatch > zMuSaMatchMapToken_
EDGetTokenT< CandidateView > zMuTkToken_
edm::Service< TFileService > fs

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 96 of file gamma_radiative_analysis.cc.

References funct::abs(), reco::Candidate::daughter(), reco::Candidate::eta(), reco::LeafCandidate::eta(), FSR_mu, FSR_mu0, FSR_mu1, FSR_tk, pat::Lepton< LeptonType >::genLepton(), pat::PATObject< ObjectType >::genParticle(), h_gamma_pt_eta_, h_mu_pt_eta_FSR_, h_mu_pt_eta_no_FSR_, i, edm::Ref< C, T, F >::isNonnull(), j, reco::Candidate::mass(), reco::Candidate::masterClone(), reco::CompositeRefCandidateT< D >::mother(), reco::Candidate::numberOfDaughters(), numOfEvent, numofGamma, reco::Candidate::pdgId(), reco::Candidate::pt(), reco::LeafCandidate::pt(), zmmcounter, zmscounter, zmtcounter, ZMuMuAnalysisNtupler_cff::zMuMu, zMuMuMatchMapToken_, zMuMuToken_, ZMuMuAnalysisNtupler_cff::zMuSa, zMuSaMatchMapToken_, zMuSaToken_, zMuTkMatchMapToken_, and zMuTkToken_.

96  {
97  Handle<CandidateView> zMuMu; //Collection of Z made by Mu global + Mu global
98  Handle<GenParticleMatch> zMuMuMatchMap; //Map of Z made by Mu global + Mu global with MC
99  event.getByToken(zMuMuToken_, zMuMu);
100  Handle<CandidateView> zMuTk; //Collection of Z made by Mu global + Track
101  Handle<GenParticleMatch> zMuTkMatchMap;
102  event.getByToken(zMuTkToken_, zMuTk);
103  Handle<CandidateView> zMuSa; //Collection of Z made by Mu global + Sa
104  Handle<GenParticleMatch> zMuSaMatchMap;
105  event.getByToken(zMuSaToken_, zMuSa);
106  numOfEvent++;
107  // ZMuMu
108  if (zMuMu->size() > 0 ) {
109  event.getByToken(zMuMuMatchMapToken_, zMuMuMatchMap);
110  for(unsigned int i = 0; i < zMuMu->size(); ++i) { //loop on candidates
111 
112  const Candidate & zMuMuCand = (*zMuMu)[i]; //the candidate
113  CandidateBaseRef zMuMuCandRef = zMuMu->refAt(i);
114 
115 
116  CandidateBaseRef dau0 = zMuMuCand.daughter(0)->masterClone();
117  CandidateBaseRef dau1 = zMuMuCand.daughter(1)->masterClone();
118  const pat::Muon& mu0 = dynamic_cast<const pat::Muon&>(*dau0);//cast in patMuon
119  const pat::Muon& mu1 = dynamic_cast<const pat::Muon&>(*dau1);
120 
121  double zmass= zMuMuCand.mass();
122  double pt0 = mu0.pt();
123  double pt1 = mu1.pt();
124  double eta0 = mu0.eta();
125  double eta1 = mu1.eta();
126  if(pt0>20 && pt1 > 20 && abs(eta0)<2 && abs(eta1)<2 && zmass > 20 && zmass < 200){
127  GenParticleRef zMuMuMatch = (*zMuMuMatchMap)[zMuMuCandRef];
128  if(zMuMuMatch.isNonnull()) { // ZMuMu matched
129  zmmcounter++;
130  FSR_mu0 = false;
131  FSR_mu1 = false;
132 
133  //MonteCarlo Study
134  const reco::GenParticle * muMc0 = mu0.genLepton();
135  const reco::GenParticle * muMc1 = mu1.genLepton();
136  const Candidate * motherMu0 = muMc0->mother();
137  const Candidate * motherMu1 = muMc1->mother();
138  int num_dau_muon0 = motherMu0->numberOfDaughters();
139  int num_dau_muon1 = motherMu1->numberOfDaughters();
140  if( num_dau_muon0 > 1 ){
141  for(int j = 0; j < num_dau_muon0; ++j){
142  int id =motherMu0 ->daughter(j)->pdgId();
143  if(id == 22){
144  double etaG = motherMu0 ->daughter(j)->eta();
145  double ptG = motherMu0 ->daughter(j)->pt();
146  h_gamma_pt_eta_->Fill(ptG,etaG);
147  h_mu_pt_eta_FSR_->Fill(pt0,eta0);
148  FSR_mu0=true;
149  numofGamma++;
150  }
151  }
152  }//end check of gamma
153  if(!FSR_mu0) h_mu_pt_eta_no_FSR_->Fill(pt0,eta0);
154  if( num_dau_muon1 > 1 ){
155  for(int j = 0; j < num_dau_muon1; ++j){
156  int id = motherMu1->daughter(j)->pdgId();
157  if(id == 22){
158  double etaG = motherMu1 ->daughter(j)->eta();
159  double ptG = motherMu1 ->daughter(j)->pt();
160  h_gamma_pt_eta_->Fill(ptG,etaG);
161  h_mu_pt_eta_FSR_->Fill(pt1,eta1);
162  FSR_mu1=true;
163  numofGamma++;
164  }
165  }
166  }//end check of gamma
167  if(!FSR_mu1) h_mu_pt_eta_no_FSR_->Fill(pt1,eta1);
168  }// end MC match
169  }//end of cuts
170  }// end loop on ZMuMu cand
171  }// end if ZMuMu size > 0
172 
173  // ZMuSa
174  if (zMuSa->size() > 0 ) {
175  event.getByToken(zMuSaMatchMapToken_, zMuSaMatchMap);
176  for(unsigned int i = 0; i < zMuSa->size(); ++i) { //loop on candidates
177 
178  const Candidate & zMuSaCand = (*zMuSa)[i]; //the candidate
179  CandidateBaseRef zMuSaCandRef = zMuSa->refAt(i);
180 
181 
182  CandidateBaseRef dau0 = zMuSaCand.daughter(0)->masterClone();
183  CandidateBaseRef dau1 = zMuSaCand.daughter(1)->masterClone();
184  const pat::Muon& mu0 = dynamic_cast<const pat::Muon&>(*dau0);//cast in patMuon
185  const pat::Muon& mu1 = dynamic_cast<const pat::Muon&>(*dau1);
186 
187  double zmass= zMuSaCand.mass();
188  double pt0 = mu0.pt();
189  double pt1 = mu1.pt();
190  double eta0 = mu0.eta();
191  double eta1 = mu1.eta();
192  if(pt0>20 && pt1 > 20 && abs(eta0)<2 && abs(eta1)<2 && zmass > 20 && zmass < 200){
193  GenParticleRef zMuSaMatch = (*zMuSaMatchMap)[zMuSaCandRef];
194  if(zMuSaMatch.isNonnull()) { // ZMuSa matched
195  FSR_mu0 = false;
196  FSR_mu1 = false;
197  zmscounter++;
198  //MonteCarlo Study
199  const reco::GenParticle * muMc0 = mu0.genLepton();
200  const reco::GenParticle * muMc1 = mu1.genLepton();
201  const Candidate * motherMu0 = muMc0->mother();
202  const Candidate * motherMu1 = muMc1->mother();
203  int num_dau_muon0 = motherMu0->numberOfDaughters();
204  int num_dau_muon1 = motherMu1->numberOfDaughters();
205  if( num_dau_muon0 > 1 ){
206  for(int j = 0; j < num_dau_muon0; ++j){
207  int id =motherMu0 ->daughter(j)->pdgId();
208  if(id == 22){
209  double etaG = motherMu0 ->daughter(j)->eta();
210  double ptG = motherMu0 ->daughter(j)->pt();
211  h_gamma_pt_eta_->Fill(ptG,etaG);
212  h_mu_pt_eta_FSR_->Fill(pt0,eta0);
213  numofGamma++;
214  FSR_mu0=true;
215  }
216  }
217  }//end check of gamma
218  if(!FSR_mu0) h_mu_pt_eta_no_FSR_->Fill(pt0,eta0);
219  if( num_dau_muon1 > 1 ){
220  for(int j = 0; j < num_dau_muon1; ++j){
221  int id = motherMu1->daughter(j)->pdgId();
222  if(id == 22){
223  double etaG = motherMu1 ->daughter(j)->eta();
224  double ptG = motherMu1 ->daughter(j)->pt();
225  h_gamma_pt_eta_->Fill(ptG,etaG);
226  h_mu_pt_eta_FSR_->Fill(pt1,eta1);
227  numofGamma++;
228  FSR_mu1=true;
229  }
230  }
231  }//end check of gamma
232  if(!FSR_mu1) h_mu_pt_eta_no_FSR_->Fill(pt1,eta1);
233  }// end MC match
234  }//end of cuts
235  }// end loop on ZMuSa cand
236  }// end if ZMuSa size > 0
237 
238 
239 
240  //ZMuTk
241  if (zMuTk->size() > 0 ) {
242  event.getByToken(zMuTkMatchMapToken_, zMuTkMatchMap);
243  for(unsigned int i = 0; i < zMuTk->size(); ++i) { //loop on candidates
244  const Candidate & zMuTkCand = (*zMuTk)[i]; //the candidate
245  CandidateBaseRef zMuTkCandRef = zMuTk->refAt(i);
246 
247 
248  CandidateBaseRef dau0 = zMuTkCand.daughter(0)->masterClone();
249  CandidateBaseRef dau1 = zMuTkCand.daughter(1)->masterClone();
250  const pat::Muon& mu0 = dynamic_cast<const pat::Muon&>(*dau0);//cast in patMuon
251  const pat::GenericParticle& mu1 = dynamic_cast<const pat::GenericParticle &>(*dau1);
252 
253 
254  double zmass= zMuTkCand.mass();
255  double pt0 = mu0.pt();
256  double pt1 = mu1.pt();
257  double eta0 = mu0.eta();
258  double eta1 = mu1.eta();
259  if(pt0>20 && pt1 > 20 && abs(eta0)<2 && abs(eta1)<2 && zmass > 20 && zmass < 200){//kinematical cuts
260  GenParticleRef zMuTkMatch = (*zMuTkMatchMap)[zMuTkCandRef];
261  if(zMuTkMatch.isNonnull()) { // ZMuTk matched
262  FSR_mu = false;
263  FSR_tk = false;
264  zmtcounter++;
265  //MonteCarlo Study
266  const reco::GenParticle * muMc0 = mu0.genLepton();
267  const reco::GenParticle * muMc1 = mu1.genParticle() ;
268  const Candidate * motherMu0 = muMc0->mother();
269  const Candidate * motherMu1 = muMc1->mother();
270  int num_dau_muon0 = motherMu0->numberOfDaughters();
271  int num_dau_muon1 = motherMu1->numberOfDaughters();
272  if( num_dau_muon0 > 1 ){
273  for(int j = 0; j < num_dau_muon0; ++j){
274  int id = motherMu0->daughter(j)->pdgId();
275  if(id == 22){
276  double etaG = motherMu0 ->daughter(j)->eta();
277  double ptG = motherMu0 ->daughter(j)->pt();
278  h_gamma_pt_eta_->Fill(ptG,etaG);
279  h_mu_pt_eta_FSR_->Fill(pt0,eta0);
280  numofGamma++;
281  FSR_mu0=true;
282  }
283  }
284  }//end check of gamma
285  if(!FSR_mu0) h_mu_pt_eta_no_FSR_->Fill(pt0,eta0);
286  if( num_dau_muon1 > 1 ){
287  for(int j = 0; j < num_dau_muon1; ++j){
288  int id = motherMu1->daughter(j)->pdgId();
289  if(id == 22){
290  double etaG = motherMu1 ->daughter(j)->eta();
291  double ptG = motherMu1 ->daughter(j)->pt();
292  h_gamma_pt_eta_->Fill(ptG,etaG);
293  h_mu_pt_eta_FSR_->Fill(pt1,eta1);
294  numofGamma++;
295  FSR_mu1=true;
296  }
297  }
298  }//end check of gamma
299  if(!FSR_mu1) h_mu_pt_eta_no_FSR_->Fill(pt1,eta1);
300  }// end MC match
301  }//end Kine-cuts
302  }// end loop on ZMuTk cand
303  }// end if ZMuTk size > 0
304 }// end analyze
int i
Definition: DBlmapReader.cc:9
EDGetTokenT< CandidateView > zMuMuToken_
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
const reco::GenParticle * genLepton() const
Definition: Lepton.h:42
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
const reco::GenParticle * genParticle(size_t idx=0) const
Definition: PATObject.h:244
virtual double pt() const =0
transverse momentum
virtual double mass() const =0
mass
Analysis-level Generic Particle class (e.g. for hadron or muon not fully reconstructed) ...
virtual size_type numberOfDaughters() const =0
number of daughters
tuple zMuSa
zMUSa vector of PSet is specific for zMuSa category
EDGetTokenT< CandidateView > zMuSaToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
EDGetTokenT< GenParticleMatch > zMuMuMatchMapToken_
EDGetTokenT< GenParticleMatch > zMuTkMatchMapToken_
virtual int pdgId() const =0
PDG identifier.
EDGetTokenT< GenParticleMatch > zMuSaMatchMapToken_
EDGetTokenT< CandidateView > zMuTkToken_
virtual double eta() const final
momentum pseudorapidity
Analysis-level muon class.
Definition: Muon.h:49
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
tuple zMuMu
zMuMu vector of PSet is common to all categories except zMuTrk category
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const final
transverse momentum
virtual const CandidateBaseRef & masterClone() const =0
void gamma_radiative_analyzer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 308 of file gamma_radiative_analysis.cc.

References gather_cfg::cout, numOfEvent, numofGamma, zmmcounter, zmscounter, and zmtcounter.

308  {
309  cout <<" ============= Summary =========="<<endl;
310  cout <<" Numero di eventi = "<< numOfEvent << endl;
311  cout <<" 1)Numero di ZMuMu matched dopo i tagli cinematici = "<< zmmcounter << endl;
312  cout <<" 2)Numero di ZMuSa matched dopo i tagli cinematici = "<< zmscounter << endl;
313  cout <<" 3)Numero di ZMuTk matched dopo i tagli cinematici = "<< zmtcounter << endl;
314  cout <<" 4)Number of gamma = "<< numofGamma << endl;
315  }
tuple cout
Definition: gather_cfg.py:145

Member Data Documentation

double gamma_radiative_analyzer::dRTrk_
private

Definition at line 65 of file gamma_radiative_analysis.cc.

double gamma_radiative_analyzer::dRVeto_
private

Definition at line 65 of file gamma_radiative_analysis.cc.

bool gamma_radiative_analyzer::FSR_mu
private

Definition at line 70 of file gamma_radiative_analysis.cc.

Referenced by analyze().

bool gamma_radiative_analyzer::FSR_mu0
private

Definition at line 70 of file gamma_radiative_analysis.cc.

Referenced by analyze().

bool gamma_radiative_analyzer::FSR_mu1
private

Definition at line 70 of file gamma_radiative_analysis.cc.

Referenced by analyze().

bool gamma_radiative_analyzer::FSR_tk
private

Definition at line 70 of file gamma_radiative_analysis.cc.

Referenced by analyze().

TH2D* gamma_radiative_analyzer::h_gamma_pt_eta_
private

Definition at line 67 of file gamma_radiative_analysis.cc.

Referenced by analyze(), and gamma_radiative_analyzer().

TH2D * gamma_radiative_analyzer::h_mu_pt_eta_FSR_
private

Definition at line 67 of file gamma_radiative_analysis.cc.

Referenced by analyze(), and gamma_radiative_analyzer().

TH2D * gamma_radiative_analyzer::h_mu_pt_eta_no_FSR_
private

Definition at line 67 of file gamma_radiative_analysis.cc.

Referenced by analyze(), and gamma_radiative_analyzer().

int gamma_radiative_analyzer::numOfEvent
private

Definition at line 73 of file gamma_radiative_analysis.cc.

Referenced by analyze(), endJob(), and gamma_radiative_analyzer().

int gamma_radiative_analyzer::numofGamma
private

Definition at line 73 of file gamma_radiative_analysis.cc.

Referenced by analyze(), endJob(), and gamma_radiative_analyzer().

double gamma_radiative_analyzer::ptThreshold_
private

Definition at line 65 of file gamma_radiative_analysis.cc.

int gamma_radiative_analyzer::zmmcounter
private

Definition at line 73 of file gamma_radiative_analysis.cc.

Referenced by analyze(), endJob(), and gamma_radiative_analyzer().

int gamma_radiative_analyzer::zmscounter
private

Definition at line 73 of file gamma_radiative_analysis.cc.

Referenced by analyze(), endJob(), and gamma_radiative_analyzer().

int gamma_radiative_analyzer::zmtcounter
private

Definition at line 73 of file gamma_radiative_analysis.cc.

Referenced by analyze(), endJob(), and gamma_radiative_analyzer().

EDGetTokenT<GenParticleMatch> gamma_radiative_analyzer::zMuMuMatchMapToken_
private

Definition at line 60 of file gamma_radiative_analysis.cc.

Referenced by analyze().

EDGetTokenT<CandidateView> gamma_radiative_analyzer::zMuMuToken_
private

Definition at line 59 of file gamma_radiative_analysis.cc.

Referenced by analyze().

EDGetTokenT<GenParticleMatch> gamma_radiative_analyzer::zMuSaMatchMapToken_
private

Definition at line 64 of file gamma_radiative_analysis.cc.

Referenced by analyze().

EDGetTokenT<CandidateView> gamma_radiative_analyzer::zMuSaToken_
private

Definition at line 63 of file gamma_radiative_analysis.cc.

Referenced by analyze().

EDGetTokenT<GenParticleMatch> gamma_radiative_analyzer::zMuTkMatchMapToken_
private

Definition at line 62 of file gamma_radiative_analysis.cc.

Referenced by analyze().

EDGetTokenT<CandidateView> gamma_radiative_analyzer::zMuTkToken_
private

Definition at line 61 of file gamma_radiative_analysis.cc.

Referenced by analyze().