CMS 3D CMS Logo

SinglePhotonJetPlusHOFilter.cc
Go to the documentation of this file.
1 // Package: SinglePhotonJetPlusHOFilter
2 // Class: SinglePhotonJetPlusHOFilter
3 //
4 /*
5  Description: [one line class summary]
6 Skimming of SinglePhoton data set for the study of HO absolute weight calculation
7 * Skimming Efficiency : ~ 2 %
8  Implementation:
9  [Notes on implementation]
10  For Secondary Datasets (SD)
11 */
12 //
13 // Original Author: Gobinda Majumder & Suman Chatterjee
14 // Created: Fri July 29 14:52:17 IST 2016
15 // $Id$
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 #include <string>
23 
24 #include <iostream>
25 #include <fstream>
26 
27 // class declaration
28 //
29 
32 
36 
40 
43 
47 
48 using namespace std;
49 using namespace edm;
50 using namespace reco;
51 
52 
54  public:
56  ~SinglePhotonJetPlusHOFilter() override;
57 
58  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
59 
60  private:
61 
62  bool filter(edm::Event&, const edm::EventSetup&) override;
63 
64  // ----------member data ---------------------------
65  int Nevt;
66  int Njetp;
67  int Npass;
68  double jtptthr;
69  double jtetath;
70  double hothres;
71  double pho_Ptcut;
73 
77 
78 };
79 
80 //
81 // constants, enums and typedefs
82 //
83 
84 //
85 // static data member definitions
86 //
87 
88 //
89 // constructors and destructor
90 //
92 
93  tok_PFJets_ = consumes<reco::PFJetCollection>( iConfig.getParameter<edm::InputTag>("PFJets"));
94  tok_hoht_ = consumes<reco::PFClusterCollection>( iConfig.getParameter<edm::InputTag>("particleFlowClusterHO"));
95  tok_photons_ = consumes<edm::View<reco::Photon> > ( iConfig.getParameter<edm::InputTag>("Photons"));
96 
97  //now do what ever initialization is needed
98  jtptthr = iConfig.getUntrackedParameter<double>("Ptcut", 90.0);
99  jtetath = iConfig.getUntrackedParameter<double>("Etacut",1.5);
100  hothres = iConfig.getUntrackedParameter<double>("HOcut", 5);
101  pho_Ptcut = iConfig.getUntrackedParameter<double>("Pho_Ptcut",120) ;
102 
103  Nevt = 0;
104  Njetp=0;
105  Npass=0;
106 }
107 
109 {
110 
111  // do anything here that needs to be done at desctruction time
112  // (e.g. close files, deallocate resources etc.)
113 
114 }
115 
116 // ------------ method called on each new Event ------------
117 bool
119 {
120  using namespace std;
121  using namespace edm;
122  using namespace reco;
123 
124  Nevt++;
125 
127  iEvent.getByToken(tok_PFJets_,PFJets);
128  bool passed=false;
129  vector <pair<double, double> > jetdirection;
130  vector<double> jetspt;
131  if (PFJets.isValid()) {
132  for (unsigned jet = 0; jet<PFJets->size(); jet++) {
133  if(((*PFJets)[jet].pt()<jtptthr)||(abs((*PFJets)[jet].eta())>jtetath)) continue;
134 
135  std::pair<double, double> etaphi((*PFJets)[jet].eta(),(*PFJets)[jet].phi());
136  jetdirection.push_back(etaphi);
137  jetspt.push_back((*PFJets)[jet].pt());
138  passed = true;
139  }
140  }
141  if (!passed) return passed;
142 
144  iEvent.getByToken(tok_photons_, photons);
145  bool pho_passed=false ;
146  vector <pair<double,double> > phodirection ;
147  vector<double> phopT ;
148 
149  if(photons.isValid()){
150 
152  for( gamma1 = photons->begin(); gamma1 != photons->end(); gamma1++ ) {
153  if(gamma1->pt()<pho_Ptcut) continue ;
154  std::pair<double,double> etaphi_pho(gamma1->eta(),gamma1->phi());
155  phodirection.push_back(etaphi_pho) ;
156  phopT.push_back(gamma1->pt()) ;
157 
158  pho_passed =true ;
159  }
160  }
161 
162  if (!pho_passed) return false;
163  Njetp++;
164  bool isJetDir=false;
165 
167  iEvent.getByToken(tok_hoht_,hoht);
168  if (hoht.isValid()) {
169  if (!(*hoht).empty()) {
170  for (unsigned ijet = 0; ijet< jetdirection.size(); ijet++) {
171 
172  bool matched=false;
173  for (unsigned iph=0; iph<phodirection.size(); iph++) {
174  if(abs(deltaPhi(phodirection[iph].second, jetdirection[ijet].second))>2.0) {
175  matched=true; break;
176  }
177  }
178  if (matched) {
179  for (PFClusterCollection::const_iterator ij=(*hoht).begin(); ij!=(*hoht).end(); ij++){
180  double hoenr = (*ij).energy();
181  if (hoenr <hothres) continue;
182 
183  const math::XYZPoint& cluster_pos = ij->position();
184 
185  double hoeta = cluster_pos.eta() ;
186  double hophi = cluster_pos.phi() ;
187 
188  double delta = deltaR2(jetdirection[ijet].first, jetdirection[ijet].second, hoeta, hophi);
189  if (delta <0.5) {
190  isJetDir=true; break;
191  }
192  }
193  }
194  if (isJetDir) { break;}
195  }
196  }
197  }
198 
199  if (isJetDir) {Npass++;}
200 
201  return isJetDir;
202 
203 }
204 
205 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
206 void
208  //The following says we do not know what parameters are allowed so do no validation
209  // Please change this to state exactly what you do use, even if it is no parameters
211  desc.setUnknown();
212  descriptions.addDefault(desc);
213 }
214 
215 //define this as a plug-in
217 
218 /*
219 
220 End of SinglePhotonJetPlusHOFilter with event 100 Jetpassed 16 passed 3
221 
222 
223 */
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
SinglePhotonJetPlusHOFilter(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::PFClusterCollection > tok_hoht_
U second(std::pair< T, U > const &p)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void addDefault(ParameterSetDescription const &psetDescription)
std::pair< T, T > etaphi(T x, T y, T z)
Definition: FastMath.h:128
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:74
bool filter(edm::Event &, const edm::EventSetup &) override
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EDGetTokenT< reco::PFJetCollection > tok_PFJets_
edm::EDGetTokenT< edm::View< reco::Photon > > tok_photons_
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86