ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitProd/TreeFiller/src/FillerGenJets.cc
Revision: 1.6
Committed: Thu Mar 18 20:21:00 2010 UTC (15 years, 1 month ago) by bendavid
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_032, Mit_031, Mit_025c_branch2, Mit_025c_branch1, Mit_030, Mit_029c, Mit_029b, Mit_030_pre1, Mit_029a, Mit_029, Mit_029_pre1, Mit_028a, Mit_025c_branch0, Mit_028, Mit_027a, Mit_027, Mit_026, Mit_025e, Mit_025d, Mit_025c, Mit_025b, Mit_025a, 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, Mit_020pre1, Mit_018, Mit_017, Mit_017pre3, Mit_017pre2, Mit_017pre1, Mit_016, Mit_015b, Mit_015a, Mit_015, Mit_014e, Mit_014d, Mit_014c, Mit_014b, Mit_014a, Mit_014, Mit_014pre3, Mit_014pre2, Mit_014pre1, Mit_013d, Mit_013c, Mit_013b, Mit_013a, Mit_013, HEAD
Branch point for: Mit_025c_branch
Changes since 1.5: +2 -2 lines
Log Message:
Fix beginrun,beginjob mess

File Contents

# User Rev Content
1 bendavid 1.6 // $Id: FillerGenJets.cc,v 1.5 2009/09/25 08:42:50 loizides Exp $
2 sixie 1.1
3     #include "MitProd/TreeFiller/interface/FillerGenJets.h"
4     #include "DataFormats/JetReco/interface/GenJet.h"
5     #include "SimDataFormats/JetMatching/interface/JetFlavour.h"
6     #include "SimDataFormats/JetMatching/interface/JetFlavourMatching.h"
7     #include "SimDataFormats/JetMatching/interface/MatchedPartons.h"
8     #include "SimDataFormats/JetMatching/interface/JetMatchedPartons.h"
9 loizides 1.4 #include "MitAna/DataTree/interface/GenJetCol.h"
10     #include "MitAna/DataTree/interface/Names.h"
11     #include "MitProd/ObjectService/interface/ObjectService.h"
12 sixie 1.1
13     using namespace std;
14     using namespace edm;
15     using namespace mithep;
16    
17     //--------------------------------------------------------------------------------------------------
18     FillerGenJets::FillerGenJets(const ParameterSet &cfg, const char *name, bool active) :
19     BaseFiller(cfg,name,active),
20     flavorMatchingActive_(Conf().getUntrackedParameter<bool>("flavorMatchingActive",true)),
21     edmName_(Conf().getUntrackedParameter<string>("edmName","genjets")),
22     mitName_(Conf().getUntrackedParameter<string>("mitName","GenJets")),
23     flavorMatchingByReferenceName_(Conf().getUntrackedParameter<string>
24     ("flavorMatchingByReferenceName","srcByReference")),
25     flavorMatchingDefinition_(Conf().getUntrackedParameter<string>
26     ("flavorMatchingDefinition","Algorithmic")),
27     genjets_(new mithep::GenJetArr(16))
28     {
29     // Constructor.
30     }
31    
32     //--------------------------------------------------------------------------------------------------
33     FillerGenJets::~FillerGenJets()
34     {
35     // Destructor.
36    
37     delete genjets_;
38     }
39    
40     //--------------------------------------------------------------------------------------------------
41 bendavid 1.6 void FillerGenJets::BookDataBlock(TreeWriter &tws)
42 sixie 1.1 {
43     // Add jets branch to tree.
44    
45 loizides 1.3 tws.AddBranch(mitName_,&genjets_);
46     OS()->add<mithep::GenJetArr>(genjets_,mitName_);
47 sixie 1.1 }
48    
49     //--------------------------------------------------------------------------------------------------
50     void FillerGenJets::FillDataBlock(const edm::Event &event,
51 loizides 1.3 const edm::EventSetup &setup)
52 sixie 1.1 {
53     // Fill jets from edm collection into our collection.
54    
55 bendavid 1.2 genjets_->Delete();
56 sixie 1.1
57 loizides 1.3 // handle for the jet collection
58 sixie 1.1 Handle<reco::GenJetCollection> hGenJetProduct;
59     GetProduct(edmName_, hGenJetProduct, event);
60    
61 loizides 1.3 // handles for jet flavour matching
62 sixie 1.1 Handle<reco::JetMatchedPartonsCollection> hPartonMatchingProduct;
63     if (flavorMatchingActive_)
64     GetProduct(flavorMatchingByReferenceName_, hPartonMatchingProduct, event);
65    
66     const reco::GenJetCollection inJets = *(hGenJetProduct.product());
67    
68     // loop through all jets
69     for (reco::GenJetCollection::const_iterator inJet = inJets.begin();
70     inJet != inJets.end(); ++inJet) {
71    
72     reco::GenJetRef jetRef(hGenJetProduct, inJet - inJets.begin());
73    
74     mithep::GenJet *jet = genjets_->Allocate();
75     new (jet) mithep::GenJet(inJet->p4().x(),
76     inJet->p4().y(),
77     inJet->p4().z(),
78     inJet->p4().e());
79    
80     jet->SetHadEnergy (inJet->hadEnergy());
81     jet->SetEmEnergy (inJet->emEnergy());
82     jet->SetInvisibleEnergy (inJet->invisibleEnergy());
83     jet->SetAuxiliaryEnergy (inJet->auxiliaryEnergy());
84    
85     if (flavorMatchingActive_) {
86     unsigned int iJet = inJet - inJets.begin();
87     const reco::JetMatchedPartonsCollection *matchedPartons = hPartonMatchingProduct.product();
88     reco::MatchedPartons matchedParton = (*matchedPartons)[matchedPartons->key(iJet)];
89     int flavorPhysDef = (matchedParton.physicsDefinitionParton().isNonnull())?
90     matchedParton.physicsDefinitionParton()->pdgId():0;
91     int flavorAlgDef = (matchedParton.algoDefinitionParton().isNonnull())?
92     matchedParton.algoDefinitionParton()->pdgId():0;
93    
94     if (flavorMatchingDefinition_ == "Algorithmic") {
95     jet->SetMatchedMCFlavor(flavorAlgDef);
96     } else if(flavorMatchingDefinition_ == "Physics") {
97     jet->SetMatchedMCFlavor(flavorPhysDef);
98     } else {
99     jet->SetMatchedMCFlavor(0);
100     }
101     }
102     }
103     genjets_->Trim();
104     }