1 |
//
|
2 |
// $Id: PATJetProducer.cc,v 1.26.2.9 2009/04/30 09:11:46 gpetrucc Exp $
|
3 |
//
|
4 |
|
5 |
#include "PhysicsTools/PatAlgos/plugins/PATJetProducer.h"
|
6 |
|
7 |
#include "FWCore/MessageLogger/interface/MessageLogger.h"
|
8 |
#include "FWCore/ParameterSet/interface/FileInPath.h"
|
9 |
|
10 |
#include "DataFormats/Common/interface/ValueMap.h"
|
11 |
#include "DataFormats/Common/interface/Association.h"
|
12 |
#include "DataFormats/Candidate/interface/CandAssociation.h"
|
13 |
|
14 |
#include "DataFormats/JetReco/interface/JetTracksAssociation.h"
|
15 |
#include "DataFormats/BTauReco/interface/JetTag.h"
|
16 |
#include "DataFormats/BTauReco/interface/TrackProbabilityTagInfo.h"
|
17 |
#include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
|
18 |
#include "DataFormats/BTauReco/interface/TrackCountingTagInfo.h"
|
19 |
#include "DataFormats/BTauReco/interface/SecondaryVertexTagInfo.h"
|
20 |
#include "DataFormats/BTauReco/interface/SoftLeptonTagInfo.h"
|
21 |
|
22 |
#include "DataFormats/Candidate/interface/CandMatchMap.h"
|
23 |
#include "SimDataFormats/JetMatching/interface/JetFlavourMatching.h"
|
24 |
|
25 |
#include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
|
26 |
#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
|
27 |
|
28 |
#include "DataFormats/Math/interface/deltaR.h"
|
29 |
|
30 |
#include "DataFormats/PatCandidates/interface/JetCorrFactors.h"
|
31 |
|
32 |
#include "FWCore/Framework/interface/Selector.h"
|
33 |
|
34 |
#include <vector>
|
35 |
#include <memory>
|
36 |
|
37 |
|
38 |
using namespace pat;
|
39 |
|
40 |
|
41 |
PATJetProducer::PATJetProducer(const edm::ParameterSet& iConfig) :
|
42 |
userDataHelper_ ( iConfig.getParameter<edm::ParameterSet>("userData") )
|
43 |
{
|
44 |
// initialize the configurables
|
45 |
jetsSrc_ = iConfig.getParameter<edm::InputTag> ( "jetSource" );
|
46 |
embedCaloTowers_ = iConfig.getParameter<bool> ( "embedCaloTowers" );
|
47 |
getJetMCFlavour_ = iConfig.getParameter<bool> ( "getJetMCFlavour" );
|
48 |
jetPartonMapSource_ = iConfig.getParameter<edm::InputTag> ( "JetPartonMapSource" );
|
49 |
addGenPartonMatch_ = iConfig.getParameter<bool> ( "addGenPartonMatch" );
|
50 |
embedGenPartonMatch_ = iConfig.getParameter<bool> ( "embedGenPartonMatch" );
|
51 |
genPartonSrc_ = iConfig.getParameter<edm::InputTag> ( "genPartonMatch" );
|
52 |
addGenJetMatch_ = iConfig.getParameter<bool> ( "addGenJetMatch" );
|
53 |
genJetSrc_ = iConfig.getParameter<edm::InputTag> ( "genJetMatch" );
|
54 |
addPartonJetMatch_ = iConfig.getParameter<bool> ( "addPartonJetMatch" );
|
55 |
partonJetSrc_ = iConfig.getParameter<edm::InputTag> ( "partonJetSource" );
|
56 |
addJetCorrFactors_ = iConfig.getParameter<bool> ( "addJetCorrFactors" );
|
57 |
jetCorrFactorsSrc_ = iConfig.getParameter<std::vector<edm::InputTag> >( "jetCorrFactorsSource" );
|
58 |
addTrigMatch_ = iConfig.getParameter<bool> ( "addTrigMatch" );
|
59 |
trigMatchSrc_ = iConfig.getParameter<std::vector<edm::InputTag> >( "trigPrimMatch" );
|
60 |
addBTagInfo_ = iConfig.getParameter<bool> ( "addBTagInfo" );
|
61 |
addDiscriminators_ = iConfig.getParameter<bool> ( "addDiscriminators" );
|
62 |
discriminatorTags_ = iConfig.getParameter<std::vector<edm::InputTag> >( "discriminatorSources" );
|
63 |
addTagInfos_ = iConfig.getParameter<bool> ( "addTagInfos" );
|
64 |
tagInfoTags_ = iConfig.getParameter<std::vector<edm::InputTag> >( "tagInfoSources" );
|
65 |
addAssociatedTracks_ = iConfig.getParameter<bool> ( "addAssociatedTracks" );
|
66 |
trackAssociation_ = iConfig.getParameter<edm::InputTag> ( "trackAssociationSource" );
|
67 |
addJetCharge_ = iConfig.getParameter<bool> ( "addJetCharge" );
|
68 |
jetCharge_ = iConfig.getParameter<edm::InputTag> ( "jetChargeSource" );
|
69 |
|
70 |
// Efficiency configurables
|
71 |
addEfficiencies_ = iConfig.getParameter<bool>("addEfficiencies");
|
72 |
if (addEfficiencies_) {
|
73 |
efficiencyLoader_ = pat::helper::EfficiencyLoader(iConfig.getParameter<edm::ParameterSet>("efficiencies"));
|
74 |
}
|
75 |
|
76 |
// Resolution configurables
|
77 |
addResolutions_ = iConfig.getParameter<bool>("addResolutions");
|
78 |
if (addResolutions_) {
|
79 |
resolutionLoader_ = pat::helper::KinResolutionsLoader(iConfig.getParameter<edm::ParameterSet>("resolutions"));
|
80 |
}
|
81 |
|
82 |
|
83 |
if (discriminatorTags_.empty()) {
|
84 |
addDiscriminators_ = false;
|
85 |
} else {
|
86 |
for (std::vector<edm::InputTag>::const_iterator it = discriminatorTags_.begin(), ed = discriminatorTags_.end(); it != ed; ++it) {
|
87 |
std::string label = it->label();
|
88 |
std::string::size_type pos = label.find("JetTags");
|
89 |
if ((pos != std::string::npos) && (pos != label.length() - 7)) {
|
90 |
label.erase(pos+7); // trim a tail after "JetTags"
|
91 |
}
|
92 |
discriminatorLabels_.push_back(label);
|
93 |
}
|
94 |
}
|
95 |
if (tagInfoTags_.empty()) {
|
96 |
addTagInfos_ = false;
|
97 |
} else {
|
98 |
for (std::vector<edm::InputTag>::const_iterator it = tagInfoTags_.begin(), ed = tagInfoTags_.end(); it != ed; ++it) {
|
99 |
std::string label = it->label();
|
100 |
std::string::size_type pos = label.find("TagInfos");
|
101 |
if ((pos != std::string::npos) && (pos != label.length() - 8)) {
|
102 |
label.erase(pos+8); // trim a tail after "TagInfos"
|
103 |
}
|
104 |
tagInfoLabels_.push_back(label);
|
105 |
}
|
106 |
}
|
107 |
|
108 |
if (!addBTagInfo_) { addDiscriminators_ = false; addTagInfos_ = false; }
|
109 |
|
110 |
// Check to see if the user wants to add user data
|
111 |
useUserData_ = false;
|
112 |
if ( iConfig.exists("userData") ) {
|
113 |
useUserData_ = true;
|
114 |
}
|
115 |
|
116 |
// produces vector of jets
|
117 |
produces<std::vector<Jet> >();
|
118 |
}
|
119 |
|
120 |
|
121 |
PATJetProducer::~PATJetProducer() {
|
122 |
|
123 |
}
|
124 |
|
125 |
|
126 |
void PATJetProducer::produce(edm::Event & iEvent, const edm::EventSetup & iSetup) {
|
127 |
|
128 |
// Get the vector of jets
|
129 |
edm::Handle<edm::View<JetType> > jets;
|
130 |
iEvent.getByLabel(jetsSrc_, jets);
|
131 |
|
132 |
if (efficiencyLoader_.enabled()) efficiencyLoader_.newEvent(iEvent);
|
133 |
if (resolutionLoader_.enabled()) resolutionLoader_.newEvent(iEvent, iSetup);
|
134 |
|
135 |
// for jet flavour
|
136 |
edm::Handle<reco::JetFlavourMatchingCollection> jetFlavMatch;
|
137 |
if (getJetMCFlavour_) iEvent.getByLabel (jetPartonMapSource_, jetFlavMatch);
|
138 |
|
139 |
// Get the vector of generated particles from the event if needed
|
140 |
edm::Handle<edm::Association<reco::GenParticleCollection> > partonMatch;
|
141 |
if (addGenPartonMatch_) iEvent.getByLabel(genPartonSrc_, partonMatch);
|
142 |
// Get the vector of GenJets from the event if needed
|
143 |
edm::Handle<edm::Association<reco::GenJetCollection> > genJetMatch;
|
144 |
if (addGenJetMatch_) iEvent.getByLabel(genJetSrc_, genJetMatch);
|
145 |
/* TO BE IMPLEMENTED FOR >= 1_5_X
|
146 |
// Get the vector of PartonJets from the event if needed
|
147 |
edm::Handle<edm::View<reco::SomePartonJetType> > partonJets;
|
148 |
if (addPartonJetMatch_) iEvent.getByLabel(partonJetSrc_, partonJets);
|
149 |
*/
|
150 |
|
151 |
// read in the jet correction factors ValueMaps
|
152 |
std::vector<edm::ValueMap<JetCorrFactors> > jetCorrs;
|
153 |
if (addJetCorrFactors_) {
|
154 |
for ( size_t i = 0; i < jetCorrFactorsSrc_.size(); ++i ) {
|
155 |
edm::Handle<edm::ValueMap<JetCorrFactors> > jetCorr;
|
156 |
iEvent.getByLabel(jetCorrFactorsSrc_[i], jetCorr);
|
157 |
jetCorrs.push_back( *jetCorr );
|
158 |
}
|
159 |
}
|
160 |
|
161 |
// Get the vector of jet tags with b-tagging info
|
162 |
std::vector<edm::Handle<reco::JetFloatAssociation::Container> > jetDiscriminators;
|
163 |
if (addBTagInfo_ && addDiscriminators_) {
|
164 |
jetDiscriminators.resize(discriminatorTags_.size());
|
165 |
for (size_t i = 0; i < discriminatorTags_.size(); ++i) {
|
166 |
iEvent.getByLabel(discriminatorTags_[i], jetDiscriminators[i]);
|
167 |
}
|
168 |
}
|
169 |
std::vector<edm::Handle<edm::View<reco::BaseTagInfo> > > jetTagInfos;
|
170 |
if (addBTagInfo_ && addTagInfos_) {
|
171 |
jetTagInfos.resize(tagInfoTags_.size());
|
172 |
for (size_t i = 0; i < tagInfoTags_.size(); ++i) {
|
173 |
iEvent.getByLabel(tagInfoTags_[i], jetTagInfos[i]);
|
174 |
}
|
175 |
}
|
176 |
|
177 |
// tracks Jet Track Association
|
178 |
edm::Handle<reco::JetTracksAssociation::Container > hTrackAss;
|
179 |
if (addAssociatedTracks_) iEvent.getByLabel(trackAssociation_, hTrackAss);
|
180 |
edm::Handle<reco::JetFloatAssociation::Container > hJetChargeAss;
|
181 |
if (addJetCharge_) iEvent.getByLabel(jetCharge_, hJetChargeAss);
|
182 |
|
183 |
// loop over jets
|
184 |
std::vector<Jet> * patJets = new std::vector<Jet>();
|
185 |
for (edm::View<JetType>::const_iterator itJet = jets->begin(); itJet != jets->end(); itJet++) {
|
186 |
|
187 |
// construct the Jet from the ref -> save ref to original object
|
188 |
unsigned int idx = itJet - jets->begin();
|
189 |
edm::RefToBase<JetType> jetRef = jets->refAt(idx);
|
190 |
edm::Ptr<JetType> jetPtr = jets->ptrAt(idx);
|
191 |
Jet ajet(jetRef);
|
192 |
|
193 |
// ensure the internal storage of the jet constituents
|
194 |
if (ajet.isCaloJet() && embedCaloTowers_) {
|
195 |
const reco::CaloJet *cj = dynamic_cast<const reco::CaloJet *>(jetRef.get());
|
196 |
ajet.setCaloTowers( cj->getCaloConstituents() );
|
197 |
}
|
198 |
|
199 |
// add jet energy scale corrections
|
200 |
if (addJetCorrFactors_) {
|
201 |
// in case only one set of jet correction factors is used, clear the string
|
202 |
// that contains the name of the jcf-module, to save storage per jet:
|
203 |
if (jetCorrFactorsSrc_.size()<=1)
|
204 |
jetCorrs.front()[jetRef].clearLabel();
|
205 |
// The default jet correction is the first in the vector
|
206 |
const JetCorrFactors & jcf = jetCorrs.front()[jetRef];
|
207 |
// uncomment for debugging
|
208 |
// jcf.print();
|
209 |
//attach first (default) jet correction factors set to the jet
|
210 |
ajet.setCorrFactors(jcf);
|
211 |
// set current default which is JetCorrFactors::L3, change P4 of ajet
|
212 |
ajet.setCorrStep(JetCorrFactors::L3);
|
213 |
|
214 |
// add additional JetCorrs for syst. studies, if present
|
215 |
for ( size_t i = 1; i < jetCorrFactorsSrc_.size(); ++i ) {
|
216 |
const JetCorrFactors & jcf = jetCorrs[i][jetRef];
|
217 |
ajet.addCorrFactors(jcf);
|
218 |
}
|
219 |
}
|
220 |
|
221 |
// get the MC flavour information for this jet
|
222 |
if (getJetMCFlavour_) {
|
223 |
ajet.setPartonFlavour( (*jetFlavMatch)[edm::RefToBase<reco::Jet>(jetRef)].getFlavour() );
|
224 |
}
|
225 |
// store the match to the generated partons
|
226 |
if (addGenPartonMatch_) {
|
227 |
reco::GenParticleRef parton = (*partonMatch)[jetRef];
|
228 |
if (parton.isNonnull() && parton.isAvailable()) {
|
229 |
ajet.setGenParton(parton, embedGenPartonMatch_);
|
230 |
} // leave empty if no match found
|
231 |
}
|
232 |
// store the match to the GenJets
|
233 |
if (addGenJetMatch_) {
|
234 |
reco::GenJetRef genjet = (*genJetMatch)[jetRef];
|
235 |
if (genjet.isNonnull() && genjet.isAvailable()) {
|
236 |
ajet.setGenJet(*genjet);
|
237 |
} // leave empty if no match found
|
238 |
}
|
239 |
|
240 |
if (efficiencyLoader_.enabled()) {
|
241 |
efficiencyLoader_.setEfficiencies( ajet, jetRef );
|
242 |
}
|
243 |
|
244 |
// IMPORTANT: DO THIS AFTER JES CORRECTIONS
|
245 |
if (resolutionLoader_.enabled()) {
|
246 |
resolutionLoader_.setResolutions(ajet);
|
247 |
}
|
248 |
|
249 |
// TO BE IMPLEMENTED FOR >=1_5_X: do the PartonJet matching
|
250 |
if (addPartonJetMatch_) {
|
251 |
}
|
252 |
|
253 |
// matches to trigger primitives
|
254 |
if ( addTrigMatch_ ) {
|
255 |
for ( size_t i = 0; i < trigMatchSrc_.size(); ++i ) {
|
256 |
edm::Handle<edm::Association<TriggerPrimitiveCollection> > trigMatch;
|
257 |
iEvent.getByLabel(trigMatchSrc_[i], trigMatch);
|
258 |
TriggerPrimitiveRef trigPrim = (*trigMatch)[jetRef];
|
259 |
if ( trigPrim.isNonnull() && trigPrim.isAvailable() ) {
|
260 |
ajet.addTriggerMatch(*trigPrim);
|
261 |
}
|
262 |
}
|
263 |
}
|
264 |
|
265 |
// add b-tag info if available & required
|
266 |
if (addBTagInfo_) {
|
267 |
if (addDiscriminators_) {
|
268 |
for (size_t k=0; k<jetDiscriminators.size(); ++k) {
|
269 |
float value = (*jetDiscriminators[k])[jetRef];
|
270 |
ajet.addBDiscriminatorPair(std::make_pair(discriminatorLabels_[k], value));
|
271 |
}
|
272 |
}
|
273 |
if (addTagInfos_) {
|
274 |
for (size_t k=0; k<jetTagInfos.size(); ++k) {
|
275 |
const edm::View<reco::BaseTagInfo> & taginfos = *jetTagInfos[k];
|
276 |
// This is not associative, so we have to search the jet
|
277 |
const reco::BaseTagInfo * match = 0;
|
278 |
// Try first by 'same index'
|
279 |
if ((idx < taginfos.size()) && (taginfos[idx].jet() == jetRef)) {
|
280 |
match = &taginfos[idx];
|
281 |
} else {
|
282 |
// otherwise fail back to a simple search
|
283 |
for (edm::View<reco::BaseTagInfo>::const_iterator itTI = taginfos.begin(), edTI = taginfos.end(); itTI != edTI; ++itTI) {
|
284 |
if (itTI->jet() == jetRef) { match = &*itTI; break; }
|
285 |
}
|
286 |
}
|
287 |
if (match != 0) ajet.addTagInfo(tagInfoLabels_[k], *match);
|
288 |
}
|
289 |
}
|
290 |
}
|
291 |
|
292 |
if (addAssociatedTracks_) ajet.setAssociatedTracks( (*hTrackAss)[jetRef] );
|
293 |
|
294 |
if (addJetCharge_) ajet.setJetCharge( (*hJetChargeAss)[jetRef] );
|
295 |
|
296 |
|
297 |
if ( useUserData_ ) {
|
298 |
userDataHelper_.add( ajet, iEvent, iSetup );
|
299 |
}
|
300 |
|
301 |
|
302 |
patJets->push_back(ajet);
|
303 |
}
|
304 |
|
305 |
// sort jets in Et
|
306 |
std::sort(patJets->begin(), patJets->end(), pTComparator_);
|
307 |
|
308 |
// put genEvt in Event
|
309 |
std::auto_ptr<std::vector<Jet> > myJets(patJets);
|
310 |
iEvent.put(myJets);
|
311 |
|
312 |
}
|
313 |
|
314 |
#include "FWCore/Framework/interface/MakerMacros.h"
|
315 |
|
316 |
DEFINE_FWK_MODULE(PATJetProducer);
|
317 |
|