CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTScoutingPFProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTrigger/JetMET
4 // Class: HLTScoutingPFProducer
5 //
11 //
12 // Original Author: Dustin James Anderson
13 // Created: Fri, 12 Jun 2015 15:49:20 GMT
14 //
15 //
16 
17 // system include files
18 #include <memory>
19 
20 // user include files
26 
35 
39 
41 
43  public:
46 
47  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
48 
49  private:
50  virtual void produce(edm::StreamID sid, edm::Event & iEvent, edm::EventSetup const & setup) const override final;
51 
58 
59  const double pfJetPtCut;
60  const double pfJetEtaCut;
61  const double pfCandidatePtCut;
62 
63  const bool doJetTags;
64  const bool doCandidates;
65  const bool doMet;
66 };
67 
68 //
69 // constructors and destructor
70 //
72  pfJetCollection_(consumes<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("pfJetCollection"))),
73  pfJetTagCollection_(consumes<reco::JetTagCollection>(iConfig.getParameter<edm::InputTag>("pfJetTagCollection"))),
74  pfCandidateCollection_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("pfCandidateCollection"))),
75  vertexCollection_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexCollection"))),
76  metCollection_(consumes<reco::METCollection>(iConfig.getParameter<edm::InputTag>("metCollection"))),
77  rho_(consumes<double>(iConfig.getParameter<edm::InputTag>("rho"))),
78  pfJetPtCut(iConfig.getParameter<double>("pfJetPtCut")),
79  pfJetEtaCut(iConfig.getParameter<double>("pfJetEtaCut")),
80  pfCandidatePtCut(iConfig.getParameter<double>("pfCandidatePtCut")),
81  doJetTags(iConfig.getParameter<bool>("doJetTags")),
82  doCandidates(iConfig.getParameter<bool>("doCandidates")),
83  doMet(iConfig.getParameter<bool>("doMet"))
84 {
85  //register products
86  produces<ScoutingPFJetCollection>();
87  produces<ScoutingParticleCollection>();
88  produces<ScoutingVertexCollection>();
89  produces<double>("rho");
90  produces<double>("pfMetPt");
91  produces<double>("pfMetPhi");
92 }
93 
95 { }
96 
97 // ------------ method called to produce the data ------------
99 {
100  using namespace edm;
101 
102  //get PF jets
105  edm::LogError ("HLTScoutingPFProducer") << "invalid collection: pfJetCollection" << "\n";
106  return;
107  }
108  //get PF jet tags
111  edm::LogError ("HLTScoutingPFProducer") << "invalid collection: pfJetTagCollection" << "\n";
112  return;
113  }
114  //get PF candidates
117  edm::LogError ("HLTScoutingPFProducer") << "invalid collection: pfCandidateCollection" << "\n";
118  return;
119  }
120  //get vertices
123  edm::LogError ("HLTScoutingPFProducer") << "invalid collection: vertexCollection" << "\n";
124  return;
125  }
126  //get rho
128  if(!iEvent.getByToken(rho_, rho)){
129  edm::LogError ("HLTScoutingPFProducer") << "invalid collection: rho" << "\n";
130  return;
131  }
132  std::auto_ptr<double> outRho(new double(*rho));
133  //get MET
135  if(doMet && !iEvent.getByToken(metCollection_, metCollection)){
136  edm::LogError ("HLTScoutingPFProducer") << "invalid collection: metCollection" << "\n";
137  return;
138  }
139 
140  //produce vertices
141  std::auto_ptr<ScoutingVertexCollection> outVertices(new ScoutingVertexCollection());
142  for(auto &vtx : *vertexCollection){
143  outVertices->emplace_back(
144  vtx.x(), vtx.y(), vtx.z(), vtx.zError()
145  );
146  }
147 
148  //produce PF candidates
149  std::auto_ptr<ScoutingParticleCollection> outPFCandidates(new ScoutingParticleCollection());
150  if(doCandidates){
151  for(auto &cand : *pfCandidateCollection){
152  if(cand.pt() > pfCandidatePtCut){
153  int vertex_index = -1;
154  int index_counter = 0;
155  double dr2 = 0.0001;
156  for (auto &vtx: *outVertices) {
157  double tmp_dr2 = pow(vtx.x() - cand.vx(), 2) + pow(vtx.y() - cand.vy(), 2)
158  + pow(vtx.z() - cand.vz(), 2);
159  if (tmp_dr2 < dr2) {
160  dr2 = tmp_dr2;
161  vertex_index = index_counter;
162  }
163  if (dr2 == 0.0)
164  break;
165  ++index_counter;
166  }
167  outPFCandidates->emplace_back(
168  cand.pt(), cand.eta(), cand.phi(), cand.mass(), cand.pdgId(), vertex_index
169  );
170  }
171  }
172  }
173 
174  //produce PF jets
175  std::auto_ptr<ScoutingPFJetCollection> outPFJets(new ScoutingPFJetCollection());
176  for(auto &jet : *pfJetCollection){
177  if(jet.pt() < pfJetPtCut || fabs(jet.eta()) > pfJetEtaCut) continue;
178  //find the jet tag corresponding to the jet
179  float tagValue = -20;
180  float minDR2 = 0.01;
181  if(doJetTags){
182  for(auto &tag : *pfJetTagCollection){
183  float dR2 = reco::deltaR2(jet, *(tag.first));
184  if(dR2 < minDR2){
185  minDR2 = dR2;
186  tagValue = tag.second;
187  }
188  }
189  }
190  //get the PF constituents of the jet
191  std::vector<int> candIndices;
192  if(doCandidates){
193  for(auto &cand : jet.getPFConstituents()){
194  if(cand->pt() > pfCandidatePtCut){
195  //search for the candidate in the collection
196  float minDR2 = 0.0001;
197  int matchIndex = -1;
198  int outIndex = 0;
199  for(auto &outCand : *outPFCandidates){
200  float dR2 = pow(cand->eta() - outCand.eta(), 2) + pow(cand->phi() - outCand.phi(), 2);
201  if(dR2 < minDR2){
202  minDR2 = dR2;
203  matchIndex = outIndex;
204  }
205  if(minDR2 == 0){
206  break;
207  }
208  outIndex++;
209  }
210  candIndices.push_back(matchIndex);
211  }
212  }
213  }
214  outPFJets->emplace_back(
215  jet.pt(), jet.eta(), jet.phi(), jet.mass(), jet.jetArea(),
216  jet.chargedHadronEnergy(), jet.neutralHadronEnergy(), jet.photonEnergy(),
217  jet.electronEnergy(), jet.muonEnergy(), jet.HFHadronEnergy(), jet.HFEMEnergy(),
218  jet.chargedHadronMultiplicity(), jet.neutralHadronMultiplicity(), jet.photonMultiplicity(),
219  jet.electronMultiplicity(), jet.muonMultiplicity(),
220  jet.HFHadronMultiplicity(), jet.HFEMMultiplicity(),
221  jet.hoEnergy(), tagValue, 0.0, candIndices
222  );
223  }
224 
225  double metPt = -999;
226  double metPhi = -999;
227  if(doMet){
228  metPt = metCollection->front().pt();
229  metPhi = metCollection->front().phi();
230  }
231  std::auto_ptr<double> outMetPt(new double(metPt));
232  std::auto_ptr<double> outMetPhi(new double(metPhi));
233 
234 
235  //put output
236  iEvent.put(outVertices);
237  iEvent.put(outPFCandidates);
238  iEvent.put(outPFJets);
239  iEvent.put(outRho, "rho");
240  iEvent.put(outMetPt, "pfMetPt");
241  iEvent.put(outMetPhi, "pfMetPhi");
242 }
243 
244 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
247  desc.add<edm::InputTag>("pfJetCollection",edm::InputTag("hltAK4PFJets"));
248  desc.add<edm::InputTag>("pfJetTagCollection",edm::InputTag("hltCombinedSecondaryVertexBJetTagsPF"));
249  desc.add<edm::InputTag>("pfCandidateCollection", edm::InputTag("hltParticleFlow"));
250  desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVertices"));
251  desc.add<edm::InputTag>("metCollection", edm::InputTag("hltPFMETProducer"));
252  desc.add<edm::InputTag>("rho", edm::InputTag("hltFixedGridRhoFastjetAll"));
253  desc.add<double>("pfJetPtCut", 20.0);
254  desc.add<double>("pfJetEtaCut", 3.0);
255  desc.add<double>("pfCandidatePtCut", 0.6);
256  desc.add<bool>("doJetTags", true);
257  desc.add<bool>("doCandidates", true);
258  desc.add<bool>("doMet", true);
259  descriptions.add("hltScoutingPFProducer", desc);
260 }
261 
262 //define this as a plug-in
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Definition: DDAxes.h:10
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
JetFloatAssociation::Container JetTagCollection
Definition: JetTag.h:18
tuple pfJetTagCollection
tuple vertexCollection
const edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidateCollection_
const edm::EDGetTokenT< reco::METCollection > metCollection_
tuple pfCandidateCollection
const edm::EDGetTokenT< reco::VertexCollection > vertexCollection_
tuple pfJetCollection
int iEvent
Definition: GenABIO.cc:230
std::vector< ScoutingVertex > ScoutingVertexCollection
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
std::vector< ScoutingParticle > ScoutingParticleCollection
virtual void produce(edm::StreamID sid, edm::Event &iEvent, edm::EventSetup const &setup) const overridefinal
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
const edm::EDGetTokenT< reco::PFJetCollection > pfJetCollection_
Collection of MET.
HLTScoutingPFProducer(const edm::ParameterSet &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const edm::EDGetTokenT< reco::JetTagCollection > pfJetTagCollection_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
tuple pfCandidatePtCut
const edm::EDGetTokenT< double > rho_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< PFJet > PFJetCollection
collection of PFJet objects
std::vector< ScoutingPFJet > ScoutingPFJetCollection
Definition: ScoutingPFJet.h:87
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40