CMS 3D CMS Logo

GenJetBCEnergyRatio.cc
Go to the documentation of this file.
1 //
2 // Plugin to store B and C ratio for a GenJet in the event
3 // Author: Attilio
4 // Date: 05.10.2007
5 //
6 
7 //=======================================================================
8 
9 // user include files
14 
21 
25 
32 
36 
37 #include <memory>
38 #include <string>
39 #include <iostream>
40 #include <vector>
41 #include <Math/VectorUtil.h>
42 #include <TMath.h>
43 
44 using namespace std;
45 using namespace reco;
46 using namespace edm;
47 using namespace ROOT::Math::VectorUtil;
48 using namespace JetMCTagUtils;
49 using namespace CandMCTagUtils;
50 
52 {
53  public:
55  ~GenJetBCEnergyRatio() override;
56 
58 
59  private:
60  void produce(StreamID, edm::Event&, const edm::EventSetup& ) const override;
62 
63 };
64 
65 //=========================================================================
66 
68 {
69  produces<JetBCEnergyRatioCollection>("bRatioCollection");
70  produces<JetBCEnergyRatioCollection>("cRatioCollection");
71  m_genjetsSrcToken = consumes< View <Jet> >(iConfig.getParameter<edm::InputTag>("genJets"));
72 }
73 
74 //=========================================================================
75 
77 {
78 }
79 
80 // ------------ method called to produce the data ------------
81 
83 {
84  Handle< View <Jet> > genjets;
85  iEvent.getByToken(m_genjetsSrcToken, genjets);
86 
88 
91 
92  if (!genjets.product()->empty()) {
93  const JetRef jj = genjets->refAt(0);
96  } else {
97  jtc1 = new JetBCEnergyRatioCollection();
98  jtc2 = new JetBCEnergyRatioCollection();
99  }
100 
101  std::unique_ptr<JetBCEnergyRatioCollection> bRatioColl(jtc1);
102  std::unique_ptr<JetBCEnergyRatioCollection> cRatioColl(jtc2);
103 
104  for( size_t j = 0; j != genjets->size(); ++j ) {
105 
106  float bRatio = EnergyRatioFromBHadrons( (*genjets)[j] );
107  float cRatio = EnergyRatioFromCHadrons( (*genjets)[j] );
108 
109  const JetRef & aJet = genjets->refAt(j) ;
110 
111  JetFloatAssociation::setValue(*bRatioColl, aJet, bRatio);
112  JetFloatAssociation::setValue(*cRatioColl, aJet, cRatio);
113 
114  }
115 
116 
117  iEvent.put(std::move(bRatioColl), "bRatioCollection");
118  iEvent.put(std::move(cRatioColl), "cRatioCollection");
119 
120 }
121 
122 //define this as a plug-in
124 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool setValue(Container &, const reco::JetBaseRef &, const JetExtendedData &)
associate jet with value. Returns false and associate nothing if jet is already associated ...
edm::Ref< JetBxCollection > JetRef
Definition: Jet.h:13
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double EnergyRatioFromBHadrons(const reco::Candidate &c)
Definition: JetMCTag.cc:12
GenJetBCEnergyRatio(const edm::ParameterSet &)
edm::EDGetTokenT< View< Jet > > m_genjetsSrcToken
T const * product() const
Definition: Handle.h:74
RefToBaseProd< T > makeRefToBaseProdFrom(RefToBase< T > const &iRef, Event const &iEvent)
fixed size matrix
HLT enums.
double EnergyRatioFromCHadrons(const reco::Candidate &c)
Definition: JetMCTag.cc:26
def move(src, dest)
Definition: eostools.py:511
void produce(StreamID, edm::Event &, const edm::EventSetup &) const override
reco::JetFloatAssociation::Container JetBCEnergyRatioCollection