ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitProd/TreeFiller/src/FillerPFJets.cc
Revision: 1.15
Committed: Tue Mar 22 19:12:34 2011 UTC (14 years, 1 month ago) by mzanetti
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_025, Mit_025pre2, Mit_024b, Mit_025pre1, Mit_024a, Mit_024, Mit_023, Mit_022a, Mit_022, Mit_020d, TMit_020d, Mit_020c, Mit_021, Mit_021pre2, Mit_021pre1, Mit_020b, Mit_020a, Mit_020
Changes since 1.14: +4 -2 lines
Log Message:
area stored by default also for uncorrected jets (will be zero)

File Contents

# User Rev Content
1 mzanetti 1.15 // $Id: FillerPFJets.cc,v 1.14 2011/03/09 17:32:53 bendavid Exp $
2 bendavid 1.1
3     #include "MitProd/TreeFiller/interface/FillerPFJets.h"
4 loizides 1.5 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
5     #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
6     #include "DataFormats/JetReco/interface/PFJet.h"
7     #include "DataFormats/BTauReco/interface/JetTag.h"
8 bendavid 1.1 #include "SimDataFormats/JetMatching/interface/JetFlavour.h"
9     #include "SimDataFormats/JetMatching/interface/JetFlavourMatching.h"
10     #include "SimDataFormats/JetMatching/interface/MatchedPartons.h"
11     #include "SimDataFormats/JetMatching/interface/JetMatchedPartons.h"
12     #include "JetMETCorrections/Objects/interface/JetCorrector.h"
13 loizides 1.5 #include "MitAna/DataTree/interface/Names.h"
14     #include "MitAna/DataTree/interface/PFJetCol.h"
15 mzanetti 1.12 #include "MitAna/DataTree/interface/PileupEnergyDensity.h"
16 loizides 1.5 #include "MitProd/ObjectService/interface/ObjectService.h"
17 bendavid 1.1
18     using namespace std;
19     using namespace edm;
20     using namespace mithep;
21    
22     //--------------------------------------------------------------------------------------------------
23     FillerPFJets::FillerPFJets(const ParameterSet &cfg, const char *name, bool active) :
24     BaseFiller(cfg,name,active),
25     flavorMatchingActive_(Conf().getUntrackedParameter<bool>("flavorMatchingActive",true)),
26     bTaggingActive_(Conf().getUntrackedParameter<bool>("bTaggingActive",true)),
27     jetToVertexActive_(Conf().getUntrackedParameter<bool>("jetToVertexActive",true)),
28     jetCorrectionsActive_(Conf().getUntrackedParameter<bool>("jetCorrectionsActive",true)),
29 mzanetti 1.13 fastJetCorrectionsActive_(Conf().getUntrackedParameter<bool>("fastJetCorrectionsActive",false)),
30 bendavid 1.1 edmName_(Conf().getUntrackedParameter<string>("edmName","recoPFJets:iterativeCone5PFJets")),
31     mitName_(Conf().getUntrackedParameter<string>("mitName","ItrCone5PFJets")),
32     jetToVertexAlphaName_(Conf().getUntrackedParameter<string>
33     ("jetToVertexAlphaName","jetToVertexAlpha")),
34     jetToVertexBetaName_(Conf().getUntrackedParameter<string>
35     ("jetToVertexBetaName","jetToVertexBetaName")),
36 mzanetti 1.13 rhoName_(Conf().getUntrackedParameter<edm::InputTag> ("rhoName")),
37 bendavid 1.1 L2JetCorrectorName_(Conf().getUntrackedParameter<string>
38     ("L2JetCorrectorName","L2JetCorrectorName")),
39     L3JetCorrectorName_(Conf().getUntrackedParameter<string>
40     ("L3JetCorrectorName","L3JetCorrectorName")),
41     flavorMatchingByReferenceName_(Conf().getUntrackedParameter<string>
42     ("flavorMatchingByReferenceName","srcByReference")),
43     flavorMatchingDefinition_(Conf().getUntrackedParameter<string>
44     ("flavorMatchingDefinition","Algorithmic")),
45     jetProbabilityBJetTagsName_(Conf().getUntrackedParameter<string>
46     ("JetProbabilityBJetTagsName","jetProbabilityBJetTags")),
47     jetBProbabilityBJetTagsName_(Conf().getUntrackedParameter<string>
48     ("JetBProbabilityBJetTagsName","jetBProbabilityBJetTags")),
49     simpleSecondaryVertexBJetTagsName_(Conf().getUntrackedParameter<string>
50     ("SimpleSecondaryVertexBJetTagsName","simpleSecondaryVertexBJetTags")),
51 bendavid 1.14 simpleSecondaryVertexHighEffBJetTagsName_(Conf().getUntrackedParameter<string>
52     ("SimpleSecondaryVertexHighEffBJetTagsName","simpleSecondaryVertexHighEffBJetTags")),
53     simpleSecondaryVertexHighPurBJetTagsName_(Conf().getUntrackedParameter<string>
54     ("SimpleSecondaryVertexHighPurBJetTagsName","simpleSecondaryVertexHighPurBJetTags")),
55 bendavid 1.1 combinedSecondaryVertexBJetTagsName_(Conf().getUntrackedParameter<string>
56     ("CombinedSecondaryVertexBJetTagsName","combinedSecondaryVertexBJetTags")),
57     combinedSecondaryVertexMVABJetTagsName_(Conf().getUntrackedParameter<string>
58     ("CombinedSecondaryVertexMVABJetTagsName","combinedSecondaryVertexMVABJetTags")),
59 bendavid 1.14 ghostTrackBJetTagsName_(Conf().getUntrackedParameter<string>
60     ("GhostTrackBJetTagsName","ghostTrackBJetTags")),
61 bendavid 1.1 trackCountingHighEffBJetTagsName_(Conf().getUntrackedParameter<string>
62     ("TrackCountingHighEffBJetTagsName","trackCountingHighEffBJetTags")),
63     trackCountingHighPurBJetTagsName_(Conf().getUntrackedParameter<string>
64     ("TrackCountingHighPurBJetTagsName","trackCountingHighPurBJetTags")),
65     softMuonBJetTagsName_(Conf().getUntrackedParameter<string>
66     ("SoftMuonBJetTagsName","softMuonBJetTags")),
67 bendavid 1.7 softMuonByIP3dBJetTagsName_(Conf().getUntrackedParameter<string>
68     ("SoftMuonByIP3dBJetTagsName","softMuonByIP3dBJetTags")),
69     softMuonByPtBJetTagsName_(Conf().getUntrackedParameter<string>
70     ("SoftMuonByPtBJetTagsName","softMuonByPtBJetTags")),
71     softElectronByIP3dBJetTagsName_(Conf().getUntrackedParameter<string>
72     ("SoftElectronByIP3dBJetTagsName","softElectronByIP3dBJetTags")),
73     softElectronByPtBJetTagsName_(Conf().getUntrackedParameter<string>
74     ("SoftElectronByPtBJetTagsName","softElectronByPtBJetTags")),
75 bendavid 1.1 pfCandMapName_(Conf().getUntrackedParameter<string>("pfCandMapName","pfCandMapName")),
76 bendavid 1.3 jetMapName_(Conf().getUntrackedParameter<string>("jetMapName","PFJetMap")),
77 bendavid 1.1 pfCandMap_(0),
78 bendavid 1.3 jetMap_(new mithep::PFJetMap),
79 bendavid 1.1 jets_(new mithep::PFJetArr(16))
80     {
81     // Constructor.
82     }
83    
84     //--------------------------------------------------------------------------------------------------
85     FillerPFJets::~FillerPFJets()
86     {
87     // Destructor.
88    
89     delete jets_;
90 bendavid 1.3 delete jetMap_;
91 bendavid 1.1 }
92    
93     //--------------------------------------------------------------------------------------------------
94 bendavid 1.10 void FillerPFJets::BookDataBlock(TreeWriter &tws)
95 bendavid 1.1 {
96     // Add jets branch to tree.
97    
98 loizides 1.2 tws.AddBranch(mitName_,&jets_);
99     OS()->add<mithep::PFJetArr>(jets_,mitName_);
100 bendavid 1.1
101 loizides 1.2 if (!pfCandMapName_.empty()) {
102     pfCandMap_ = OS()->get<PFCandidateMap>(pfCandMapName_);
103     if (pfCandMap_)
104     AddBranchDep(mitName_,pfCandMap_->GetBrName());
105     }
106 bendavid 1.3 if (!jetMapName_.empty()) {
107     jetMap_->SetBrName(mitName_);
108     OS()->add<PFJetMap>(jetMap_,jetMapName_);
109     }
110 bendavid 1.1 }
111    
112     //--------------------------------------------------------------------------------------------------
113     void FillerPFJets::FillDataBlock(const edm::Event &event,
114 loizides 1.2 const edm::EventSetup &setup)
115 bendavid 1.1 {
116     // Fill jets from edm collection into our collection.
117    
118     jets_->Delete();
119 bendavid 1.3 jetMap_->Reset();
120 bendavid 1.1
121 loizides 1.2 // handle for the Jet Collection
122 bendavid 1.1 Handle<reco::PFJetCollection> hJetProduct;
123     GetProduct(edmName_, hJetProduct, event);
124    
125 mzanetti 1.13 // Handle<double> rho;
126     // GetProduct(rhoName_, rho, event);
127    
128 mzanetti 1.12 Handle<double> rho;
129 mzanetti 1.13 event.getByLabel(rhoName_,rho);
130 mzanetti 1.12
131 loizides 1.2 // handles for jet flavour matching
132 bendavid 1.1 Handle<reco::JetMatchedPartonsCollection> hPartonMatchingProduct;
133     if (flavorMatchingActive_)
134     GetProduct(flavorMatchingByReferenceName_, hPartonMatchingProduct, event);
135    
136     Handle<reco::JetTagCollection> hJetProbabilityBJetTags;
137     Handle<reco::JetTagCollection> hJetBProbabilityBJetTags;
138     Handle<reco::JetTagCollection> hSimpleSecondaryVertexBJetTags;
139 bendavid 1.14 Handle<reco::JetTagCollection> hSimpleSecondaryVertexHighEffBJetTags;
140     Handle<reco::JetTagCollection> hSimpleSecondaryVertexHighPurBJetTags;
141 bendavid 1.1 Handle<reco::JetTagCollection> hCombinedSecondaryVertexBJetTags;
142     Handle<reco::JetTagCollection> hCombinedSecondaryVertexMVABJetTags;
143     Handle<reco::JetTagCollection> hTrackCountingHighEffBJetTags;
144     Handle<reco::JetTagCollection> hTrackCountingHighPurBJetTags;
145     Handle<reco::JetTagCollection> hSoftMuonBJetTags;
146 bendavid 1.7 Handle<reco::JetTagCollection> hSoftMuonByIP3dBJetTags;
147     Handle<reco::JetTagCollection> hSoftMuonByPtBJetTags;
148     Handle<reco::JetTagCollection> hSoftElectronByIP3dBJetTags;
149     Handle<reco::JetTagCollection> hSoftElectronByPtBJetTags;
150 bendavid 1.14 Handle<reco::JetTagCollection> hGhostTrackBJetTags;
151 bendavid 1.1
152     if (bTaggingActive_) {
153     GetProduct(jetProbabilityBJetTagsName_, hJetProbabilityBJetTags, event);
154 bendavid 1.14 GetProduct(jetBProbabilityBJetTagsName_, hJetBProbabilityBJetTags, event);
155     event.getByLabel(simpleSecondaryVertexBJetTagsName_,hSimpleSecondaryVertexBJetTags);
156     event.getByLabel(simpleSecondaryVertexHighEffBJetTagsName_,hSimpleSecondaryVertexHighEffBJetTags);
157     event.getByLabel(simpleSecondaryVertexHighPurBJetTagsName_,hSimpleSecondaryVertexHighPurBJetTags);
158 bendavid 1.1 GetProduct(combinedSecondaryVertexBJetTagsName_, hCombinedSecondaryVertexBJetTags, event);
159 loizides 1.8 GetProduct(combinedSecondaryVertexMVABJetTagsName_, hCombinedSecondaryVertexMVABJetTags, event);
160 bendavid 1.1 GetProduct(trackCountingHighEffBJetTagsName_, hTrackCountingHighEffBJetTags, event);
161     GetProduct(trackCountingHighPurBJetTagsName_, hTrackCountingHighPurBJetTags, event);
162     GetProduct(softMuonBJetTagsName_, hSoftMuonBJetTags, event);
163 bendavid 1.7 GetProduct(softMuonByIP3dBJetTagsName_, hSoftMuonByIP3dBJetTags, event);
164     GetProduct(softMuonByPtBJetTagsName_, hSoftMuonByPtBJetTags, event);
165     GetProduct(softElectronByIP3dBJetTagsName_, hSoftElectronByIP3dBJetTags, event);
166     GetProduct(softElectronByPtBJetTagsName_, hSoftElectronByPtBJetTags, event);
167 bendavid 1.14 event.getByLabel(ghostTrackBJetTagsName_,hGhostTrackBJetTags);
168 bendavid 1.1 }
169    
170     const reco::PFJetCollection inJets = *(hJetProduct.product());
171    
172     //Handles to Jet to Vertex Association
173     Handle<std::vector<double> > JV_alpha;
174     Handle<std::vector<double> > JV_beta;
175     std::vector<double>::const_iterator it_jv_alpha;
176     std::vector<double>::const_iterator it_jv_beta;
177    
178     if (jetToVertexActive_) {
179     GetProduct(jetToVertexAlphaName_, JV_alpha, event);
180     GetProduct(jetToVertexBetaName_, JV_beta, event);
181     it_jv_alpha = JV_alpha->begin();
182     it_jv_beta = JV_beta->begin();
183     }
184    
185     //Define Jet Correction Services
186     const JetCorrector* correctorL2 = 0;
187     const JetCorrector* correctorL3 = 0;
188     if (jetCorrectionsActive_) {
189     correctorL2 = JetCorrector::getJetCorrector (L2JetCorrectorName_,setup);
190     correctorL3 = JetCorrector::getJetCorrector (L3JetCorrectorName_,setup);
191     }
192    
193     // loop through all jets
194     for (reco::PFJetCollection::const_iterator inJet = inJets.begin();
195     inJet != inJets.end(); ++inJet) {
196    
197     reco::PFJetRef jetRef(hJetProduct, inJet - inJets.begin());
198     reco::JetBaseRef jetBaseRef(jetRef);
199    
200     mithep::PFJet *jet = jets_->Allocate();
201     new (jet) mithep::PFJet(inJet->p4().x(),
202     inJet->p4().y(),
203     inJet->p4().z(),
204     inJet->p4().e());
205    
206 bendavid 1.3 //add to map
207     edm::Ptr<reco::PFJet> thePtr(hJetProduct, inJet - inJets.begin());
208     jetMap_->Add(thePtr, jet);
209 sixie 1.11
210     //fill jet moments
211     jet->SetSigmaEta(TMath::Sqrt(inJet->etaetaMoment()));
212     jet->SetSigmaPhi(TMath::Sqrt(inJet->phiphiMoment()));
213    
214    
215 bendavid 1.1 //fill pfjet-specific quantities
216     jet->SetChargedHadronEnergy(inJet->chargedHadronEnergy());
217     jet->SetNeutralHadronEnergy(inJet->neutralHadronEnergy());
218     jet->SetChargedEmEnergy(inJet->chargedEmEnergy());
219     jet->SetChargedMuEnergy(inJet->chargedMuEnergy());
220     jet->SetNeutralEmEnergy(inJet->neutralEmEnergy());
221     jet->SetChargedMultiplicity(inJet->chargedMultiplicity());
222     jet->SetNeutralMultiplicity(inJet->neutralMultiplicity());
223     jet->SetMuonMultiplicity(inJet->muonMultiplicity());
224    
225     if (jetToVertexActive_) {
226     //compute alpha and beta parameter for jets
227     jet->SetAlpha((*it_jv_alpha));
228     jet->SetBeta((*it_jv_beta));
229     }
230    
231 mzanetti 1.15 // fill the area anyway
232     jet->SetJetArea(inJet->jetArea());
233    
234 bendavid 1.1 //Jet Corrections
235 mzanetti 1.12 if (fastJetCorrectionsActive_) {
236 mzanetti 1.13 double l1Scale = (inJet->pt() - (*rho)*inJet->jetArea())/inJet->pt();
237     l1Scale = (l1Scale>0) ? l1Scale : 0.0;
238     jet->SetL1OffsetCorrectionScale( l1Scale);
239 mzanetti 1.12 jet->EnableCorrection(mithep::PFJet::L1);
240     }
241 bendavid 1.1 if (jetCorrectionsActive_) {
242     double L2Scale = correctorL2->correction(inJet->p4());
243     double L3Scale = correctorL3->correction(inJet->p4()*L2Scale);
244     jet->SetL2RelativeCorrectionScale(L2Scale);
245     jet->SetL3AbsoluteCorrectionScale(L3Scale);
246     jet->EnableCorrection(mithep::PFJet::L2);
247     jet->EnableCorrection(mithep::PFJet::L3);
248     }
249    
250     if (bTaggingActive_) {
251     jet->SetJetProbabilityBJetTagsDisc((*(hJetProbabilityBJetTags.product()))[jetBaseRef]);
252     jet->SetJetBProbabilityBJetTagsDisc((*(hJetBProbabilityBJetTags.product()))[jetBaseRef]);
253 bendavid 1.14 // jet->SetSimpleSecondaryVertexBJetTagsDisc((*(hSimpleSecondaryVertexBJetTags.product()))[jetBaseRef]);
254     jet->SetSimpleSecondaryVertexHighEffBJetTagsDisc((*(hSimpleSecondaryVertexHighEffBJetTags.product()))[jetBaseRef]);
255     jet->SetSimpleSecondaryVertexHighPurBJetTagsDisc((*(hSimpleSecondaryVertexHighPurBJetTags.product()))[jetBaseRef]);
256 loizides 1.2 jet->SetCombinedSecondaryVertexBJetTagsDisc(
257     (*(hCombinedSecondaryVertexBJetTags.product()))[jetBaseRef]);
258     jet->SetCombinedSecondaryVertexMVABJetTagsDisc(
259     (*(hCombinedSecondaryVertexMVABJetTags.product()))[jetBaseRef]);
260     jet->SetTrackCountingHighEffBJetTagsDisc(
261     (*(hTrackCountingHighEffBJetTags.product()))[jetBaseRef]);
262     jet->SetTrackCountingHighPurBJetTagsDisc(
263     (*(hTrackCountingHighPurBJetTags.product()))[jetBaseRef]);
264 bendavid 1.1 jet->SetSoftMuonBJetTagsDisc((*(hSoftMuonBJetTags.product()))[jetBaseRef]);
265 bendavid 1.7 jet->SetSoftMuonByIP3dBJetTagsDisc((*(hSoftMuonByIP3dBJetTags.product()))[jetBaseRef]);
266     jet->SetSoftMuonByPtBJetTagsDisc((*(hSoftMuonByPtBJetTags.product()))[jetBaseRef]);
267 loizides 1.8 jet->SetSoftElectronByIP3dBJetTagsDisc((*(hSoftElectronByIP3dBJetTags.product()))[jetBaseRef]);
268 bendavid 1.7 jet->SetSoftElectronByPtBJetTagsDisc((*(hSoftElectronByPtBJetTags.product()))[jetBaseRef]);
269 bendavid 1.14 // jet->SetGhostTrackBJetTagsDisc((*(hGhostTrackBJetTags.product()))[jetBaseRef]);
270 bendavid 1.1 }
271    
272 loizides 1.4 // get the Monte Carlo flavour matching information
273 bendavid 1.1 if (flavorMatchingActive_) {
274     unsigned int iJet = inJet - inJets.begin();
275     const reco::JetMatchedPartonsCollection *matchedPartons = hPartonMatchingProduct.product();
276     reco::MatchedPartons matchedParton = (*matchedPartons)[matchedPartons->key(iJet)];
277     int flavorPhysDef = (matchedParton.physicsDefinitionParton().isNonnull())?
278     matchedParton.physicsDefinitionParton()->pdgId():0;
279     int flavorAlgDef = (matchedParton.algoDefinitionParton().isNonnull())?
280     matchedParton.algoDefinitionParton()->pdgId():0;
281    
282     if (flavorMatchingDefinition_ == "Algorithmic") {
283     jet->SetMatchedMCFlavor(flavorAlgDef);
284     } else if(flavorMatchingDefinition_ == "Physics") {
285     jet->SetMatchedMCFlavor(flavorPhysDef);
286     } else {
287     jet->SetMatchedMCFlavor(0);
288     }
289     }
290    
291 loizides 1.4 // add PFCandidate refs
292 bendavid 1.1 if (pfCandMap_) {
293     for (uint i=0; i<inJet->numberOfDaughters(); ++i) {
294     const reco::CandidatePtr candPtr = inJet->daughterPtr(i);
295     const reco::PFCandidatePtr pfPtr(candPtr);
296     const mithep::PFCandidate *constituent = pfCandMap_->GetMit(pfPtr);
297     jet->AddPFCand(constituent);
298     }
299     }
300     if (jetToVertexActive_) {
301     it_jv_alpha++;
302     it_jv_beta++;
303     }
304     }
305     jets_->Trim();
306     }