CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoSelector.cc
Go to the documentation of this file.
1 /* \class RecoSelector
2 *
3 * Class to apply analysis cuts in the TriggerValidation Code
4 *
5 * Author: Massimiliano Chiorboli Date: August 2007
6 // Maurizio Pierini
7 // Maria Spiropulu
8 *
9 */
10 
13 
15 
16 using namespace edm;
17 using namespace reco;
18 using namespace std;
19 
21 {
22 
23  name = userCut_params.getParameter<string>("name");
24  m_electronSrc = userCut_params.getParameter<string>("electrons");
25  m_muonSrc = userCut_params.getParameter<string>("muons");
26  m_jetsSrc = userCut_params.getParameter<string>("jets");
27  m_photonProducerSrc = userCut_params.getParameter<string>("photonProducer");
28  m_photonSrc = userCut_params.getParameter<string>("photons");
29  m_calometSrc = userCut_params.getParameter<string>("calomet");
30 
31  reco_metMin = userCut_params.getParameter<double>("reco_metMin") ;
32  reco_ptJet1Min = userCut_params.getParameter<double>("reco_ptJet1Min");
33  reco_ptJet2Min = userCut_params.getParameter<double>("reco_ptJet2Min");
34  reco_ptElecMin = userCut_params.getParameter<double>("reco_ptElecMin");
35  reco_ptMuonMin = userCut_params.getParameter<double>("reco_ptMuonMin");
36  reco_ptPhotMin = userCut_params.getParameter<double>("reco_ptPhotMin");
37 
38  edm::LogInfo("RecoSelectorParameters") << endl;
39  edm::LogInfo("RecoSelectorParameters") << "UserAnalysis parameters for " << name << " selection:" << endl;
40  edm::LogInfo("RecoSelectorParameters") << " reco_metMin = " << reco_metMin << endl;
41  edm::LogInfo("RecoSelectorParameters") << " reco_ptJet1Min = " << reco_ptJet1Min << endl;
42  edm::LogInfo("RecoSelectorParameters") << " reco_ptJet2Min = " << reco_ptJet2Min << endl;
43  edm::LogInfo("RecoSelectorParameters") << " reco_ptElecMin = " << reco_ptElecMin << endl;
44  edm::LogInfo("RecoSelectorParameters") << " reco_ptMuonMin = " << reco_ptMuonMin << endl;
45  edm::LogInfo("RecoSelectorParameters") << " reco_ptPhotMin = " << reco_ptPhotMin << endl;
46 
47 }
48 
49 string RecoSelector::GetName(){return name;}
50 
52 {
53 
54 
55  this->handleObjects(iEvent);
56 
57 
58  bool TotalCutPassed = false;
59 
60  bool ElectronCutPassed = false;
61  for(unsigned int i=0; i<theElectronCollection->size(); i++) {
62  GsfElectron electron = (*theElectronCollection)[i];
63  if (electron.ecalDrivenSeed()) {
64  float elenergy = electron.superCluster()->energy();
65  float elpt = electron.pt() * elenergy / electron.p();
66  if(elpt>reco_ptElecMin) ElectronCutPassed = true;
67  LogDebug("RecoSelectorCuts") << "elpt = " << elpt << endl;
68  }
69  else {
70  // float elenergy = 0;
71  // float elpt = 0;
72  ElectronCutPassed = false;
73  }
74  }
75 
76  bool MuonCutPassed = false;
77  for(unsigned int i=0; i<theMuonCollection->size(); i++) {
78  Muon muon = (*theMuonCollection)[i];
79  float muonpt = muon.pt();
80  if(muonpt>reco_ptMuonMin) MuonCutPassed = true;
81  LogDebug("RecoSelectorCuts") << "muonpt = " << muonpt << endl;
82  }
83 
84  bool PhotonCutPassed = false;
85  for(unsigned int i=0; i<thePhotonCollection->size(); i++) {
86  Photon photon = (*thePhotonCollection)[i];
87  float photonpt = photon.pt();
88  if(photonpt>reco_ptPhotMin) PhotonCutPassed = true;
89  LogDebug("RecoSelectorCuts") << "photonpt = " << photonpt << endl;
90  }
91 
92 
93  bool JetCutPassed = false;
94  if(theCaloJetCollection->size()>1) {
95  if((*theCaloJetCollection)[0].pt() > reco_ptJet1Min &&
96  (*theCaloJetCollection)[1].pt() > reco_ptJet2Min ) JetCutPassed = true;
97  }
98 
99  bool MetCutPassed = false;
100  const CaloMET theCaloMET = theCaloMETCollection->front();
101  float calomet = theCaloMET.pt();
102  LogDebug("RecoSelectorCuts") << "calomet = " << calomet << endl;
103  if(calomet > reco_metMin) MetCutPassed = true;
104 
105  edm::LogInfo("RecoSelectorPassed") << "ElectronCutPassed = " << (int)ElectronCutPassed << endl;
106  edm::LogInfo("RecoSelectorPassed") << "MuonCutPassed = " << (int)MuonCutPassed << endl;
107  edm::LogInfo("RecoSelectorPassed") << "PhotonCutPassed = " << (int)PhotonCutPassed << endl;
108  edm::LogInfo("RecoSelectorPassed") << "JetCutPassed = " << (int)JetCutPassed << endl;
109  edm::LogInfo("RecoSelectorPassed") << "MetCutPassed = " << (int)MetCutPassed << endl;
110 
111 
112  if(
113  (ElectronCutPassed || MuonCutPassed) &&
114  PhotonCutPassed &&
115  JetCutPassed &&
116  MetCutPassed ) TotalCutPassed = true;
117 
118  return TotalCutPassed;
119 }
120 
122 {
123 
124  //Get the electrons
125  Handle<GsfElectronCollection> theElectronCollectionHandle;
126  iEvent.getByLabel(m_electronSrc, theElectronCollectionHandle);
127  theElectronCollection = theElectronCollectionHandle.product();
128 
129  //Get the Muons
130  Handle<MuonCollection> theMuonCollectionHandle;
131  iEvent.getByLabel(m_muonSrc, theMuonCollectionHandle);
132  theMuonCollection = theMuonCollectionHandle.product();
133 
134  //Get the Photons
135  Handle<PhotonCollection> thePhotonCollectionHandle;
136  iEvent.getByLabel(m_photonProducerSrc, m_photonSrc, thePhotonCollectionHandle);
137  thePhotonCollection = thePhotonCollectionHandle.product();
138 
139  //Get the CaloJets
140  Handle<CaloJetCollection> theCaloJetCollectionHandle;
141  iEvent.getByLabel(m_jetsSrc, theCaloJetCollectionHandle);
142  theCaloJetCollection = theCaloJetCollectionHandle.product();
143 
144  //Get the CaloMET
145  Handle<CaloMETCollection> theCaloMETCollectionHandle;
146  iEvent.getByLabel(m_calometSrc, theCaloMETCollectionHandle);
147  theCaloMETCollection = theCaloMETCollectionHandle.product();
148 
149 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
RecoSelector(edm::ParameterSet userCut_params)
Definition: RecoSelector.cc:20
virtual double p() const
magnitude of momentum vector
std::string GetName()
Definition: RecoSelector.cc:49
void handleObjects(const edm::Event &)
int iEvent
Definition: GenABIO.cc:243
virtual SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:168
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
bool isSelected(const edm::Event &)
Definition: RecoSelector.cc:51
virtual double pt() const
transverse momentum
T const * product() const
Definition: Handle.h:74
bool ecalDrivenSeed() const
Definition: GsfElectron.h:172