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 = iC.consumes<GsfElectronCollection>(userCut_params.getParameter<string>("electrons"));
25  m_muonSrc = iC.consumes<MuonCollection>(userCut_params.getParameter<string>("muons"));
26  m_jetsSrc = iC.consumes<CaloJetCollection>(userCut_params.getParameter<string>("jets"));
27  m_photonSrc = userCut_params.getParameter<string>("photons");
28  m_photonProducerSrc = userCut_params.getParameter<string>("photonProducer");
29  m_photon_token_ = iC.consumes<PhotonCollection>(edm::InputTag(m_photonProducerSrc, m_photonSrc));
30  m_calometSrc = iC.consumes<CaloMETCollection>(userCut_params.getParameter<string>("calomet"));
31 
32  reco_metMin = userCut_params.getParameter<double>("reco_metMin") ;
33  reco_ptJet1Min = userCut_params.getParameter<double>("reco_ptJet1Min");
34  reco_ptJet2Min = userCut_params.getParameter<double>("reco_ptJet2Min");
35  reco_ptElecMin = userCut_params.getParameter<double>("reco_ptElecMin");
36  reco_ptMuonMin = userCut_params.getParameter<double>("reco_ptMuonMin");
37  reco_ptPhotMin = userCut_params.getParameter<double>("reco_ptPhotMin");
38 
39  edm::LogInfo("RecoSelectorParameters") << endl;
40  edm::LogInfo("RecoSelectorParameters") << "UserAnalysis parameters for " << name << " selection:" << endl;
41  edm::LogInfo("RecoSelectorParameters") << " reco_metMin = " << reco_metMin << endl;
42  edm::LogInfo("RecoSelectorParameters") << " reco_ptJet1Min = " << reco_ptJet1Min << endl;
43  edm::LogInfo("RecoSelectorParameters") << " reco_ptJet2Min = " << reco_ptJet2Min << endl;
44  edm::LogInfo("RecoSelectorParameters") << " reco_ptElecMin = " << reco_ptElecMin << endl;
45  edm::LogInfo("RecoSelectorParameters") << " reco_ptMuonMin = " << reco_ptMuonMin << endl;
46  edm::LogInfo("RecoSelectorParameters") << " reco_ptPhotMin = " << reco_ptPhotMin << endl;
47 
48 }
49 
50 string RecoSelector::GetName(){return name;}
51 
53 {
54 
55 
56  this->handleObjects(iEvent);
57 
58 
59  bool TotalCutPassed = false;
60 
61  bool ElectronCutPassed = false;
62  for(unsigned int i=0; i<theElectronCollection->size(); i++) {
63  GsfElectron electron = (*theElectronCollection)[i];
64  if (electron.ecalDrivenSeed()) {
65  float elenergy = electron.superCluster()->energy();
66  float elpt = electron.pt() * elenergy / electron.p();
67  if(elpt>reco_ptElecMin) ElectronCutPassed = true;
68  LogDebug("RecoSelectorCuts") << "elpt = " << elpt << endl;
69  }
70  else {
71  // float elenergy = 0;
72  // float elpt = 0;
73  ElectronCutPassed = false;
74  }
75  }
76 
77  bool MuonCutPassed = false;
78  for(unsigned int i=0; i<theMuonCollection->size(); i++) {
79  Muon muon = (*theMuonCollection)[i];
80  float muonpt = muon.pt();
81  if(muonpt>reco_ptMuonMin) MuonCutPassed = true;
82  LogDebug("RecoSelectorCuts") << "muonpt = " << muonpt << endl;
83  }
84 
85  bool PhotonCutPassed = false;
86  for(unsigned int i=0; i<thePhotonCollection->size(); i++) {
87  Photon photon = (*thePhotonCollection)[i];
88  float photonpt = photon.pt();
89  if(photonpt>reco_ptPhotMin) PhotonCutPassed = true;
90  LogDebug("RecoSelectorCuts") << "photonpt = " << photonpt << endl;
91  }
92 
93 
94  bool JetCutPassed = false;
95  if(theCaloJetCollection->size()>1) {
96  if((*theCaloJetCollection)[0].pt() > reco_ptJet1Min &&
97  (*theCaloJetCollection)[1].pt() > reco_ptJet2Min ) JetCutPassed = true;
98  }
99 
100  bool MetCutPassed = false;
101  const CaloMET theCaloMET = theCaloMETCollection->front();
102  float calomet = theCaloMET.pt();
103  LogDebug("RecoSelectorCuts") << "calomet = " << calomet << endl;
104  if(calomet > reco_metMin) MetCutPassed = true;
105 
106  edm::LogInfo("RecoSelectorPassed") << "ElectronCutPassed = " << (int)ElectronCutPassed << endl;
107  edm::LogInfo("RecoSelectorPassed") << "MuonCutPassed = " << (int)MuonCutPassed << endl;
108  edm::LogInfo("RecoSelectorPassed") << "PhotonCutPassed = " << (int)PhotonCutPassed << endl;
109  edm::LogInfo("RecoSelectorPassed") << "JetCutPassed = " << (int)JetCutPassed << endl;
110  edm::LogInfo("RecoSelectorPassed") << "MetCutPassed = " << (int)MetCutPassed << endl;
111 
112 
113  if(
114  (ElectronCutPassed || MuonCutPassed) &&
115  PhotonCutPassed &&
116  JetCutPassed &&
117  MetCutPassed ) TotalCutPassed = true;
118 
119  return TotalCutPassed;
120 }
121 
123 {
124 
125  //Get the electrons
126  Handle<GsfElectronCollection> theElectronCollectionHandle;
127  iEvent.getByToken(m_electronSrc, theElectronCollectionHandle);
128  theElectronCollection = theElectronCollectionHandle.product();
129 
130  //Get the Muons
131  Handle<MuonCollection> theMuonCollectionHandle;
132  iEvent.getByToken(m_muonSrc, theMuonCollectionHandle);
133  theMuonCollection = theMuonCollectionHandle.product();
134 
135  //Get the Photons
136  Handle<PhotonCollection> thePhotonCollectionHandle;
137  iEvent.getByToken(m_photon_token_, thePhotonCollectionHandle);
138  thePhotonCollection = thePhotonCollectionHandle.product();
139 
140  //Get the CaloJets
141  Handle<CaloJetCollection> theCaloJetCollectionHandle;
142  iEvent.getByToken(m_jetsSrc, theCaloJetCollectionHandle);
143  theCaloJetCollection = theCaloJetCollectionHandle.product();
144 
145  //Get the CaloMET
146  Handle<CaloMETCollection> theCaloMETCollectionHandle;
147  iEvent.getByToken(m_calometSrc, theCaloMETCollectionHandle);
148  theCaloMETCollection = theCaloMETCollectionHandle.product();
149 
150 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual double p() const
magnitude of momentum vector
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
std::string GetName()
Definition: RecoSelector.cc:50
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
virtual double pt() const
transverse momentum
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
void handleObjects(const edm::Event &)
Collection of Calo MET.
int iEvent
Definition: GenABIO.cc:230
virtual SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:182
bool isSelected(const edm::Event &)
Definition: RecoSelector.cc:52
T const * product() const
Definition: Handle.h:81
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
RecoSelector(const edm::ParameterSet &userCut_params, edm::ConsumesCollector &&iC)
Definition: RecoSelector.cc:20
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
bool ecalDrivenSeed() const
Definition: GsfElectron.h:186