CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
L1GenTreeProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1TriggerDPG/L1Ntuples
4 // Class: L1GenTreeProducer
5 //
13 //
14 // Original Author:
15 // Created:
16 // $Id: L1GenTreeProducer.cc,v 1.8 2012/08/29 12:44:03 jbrooke Exp $
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 
23 // framework
30 
31 // input data formats
34 
37 #include "HepMC/GenParticle.h"
38 #include "HepMC/GenVertex.h"
40 
42 
43 // ROOT output stuff
46 #include "TTree.h"
47 
49 
50 //
51 // class declaration
52 //
53 
55 public:
56  explicit L1GenTreeProducer(const edm::ParameterSet&);
57  ~L1GenTreeProducer() override;
58 
59 private:
60  void beginJob(void) override;
61  void analyze(const edm::Event&, const edm::EventSetup&) override;
62  void endJob() override;
63 
64 private:
65  unsigned maxL1Upgrade_;
66 
67  // output file
69 
70  // tree
71  TTree* tree_;
72 
73  // data format
74  std::unique_ptr<L1Analysis::L1AnalysisGeneratorDataFormat> l1GenData_;
75 
76  // EDM input tags
81 };
82 
84  genJetToken_ = consumes<reco::GenJetCollection>(iConfig.getUntrackedParameter<edm::InputTag>("genJetToken"));
86  consumes<reco::GenParticleCollection>(iConfig.getUntrackedParameter<edm::InputTag>("genParticleToken"));
88  consumes<std::vector<PileupSummaryInfo>>(iConfig.getUntrackedParameter<edm::InputTag>("pileupInfoToken"));
89  genInfoToken_ = consumes<GenEventInfoProduct>(iConfig.getParameter<edm::InputTag>("genInfoToken"));
90 
91  l1GenData_ = std::make_unique<L1Analysis::L1AnalysisGeneratorDataFormat>();
92 
93  // set up output
94  tree_ = fs_->make<TTree>("L1GenTree", "L1GenTree");
95  tree_->Branch("Generator", "L1Analysis::L1AnalysisGeneratorDataFormat", l1GenData_.get(), 32000, 3);
96 }
97 
99  // do anything here that needs to be done at desctruction time
100  // (e.g. close files, deallocate resources etc.)
101 }
102 
103 //
104 // member functions
105 //
106 
107 // ------------ method called to for each event ------------
110  iEvent.getByToken(genInfoToken_, genInfo);
111 
112  l1GenData_->Reset();
113 
114  if (genInfo.isValid()) {
115  l1GenData_->weight = genInfo->weight();
116  l1GenData_->pthat = genInfo->hasBinningValues() ? (genInfo->binningValues())[0] : 0.0;
117  }
118 
120  iEvent.getByToken(genJetToken_, genJets);
121 
122  if (genJets.isValid()) {
123  reco::GenJetCollection::const_iterator jetItr = genJets->begin();
124  reco::GenJetCollection::const_iterator jetEnd = genJets->end();
125  for (; jetItr != jetEnd; ++jetItr) {
126  l1GenData_->jetPt.push_back(jetItr->pt());
127  l1GenData_->jetEta.push_back(jetItr->eta());
128  l1GenData_->jetPhi.push_back(jetItr->phi());
129  l1GenData_->nJet++;
130  }
131 
132  } else {
133  edm::LogWarning("MissingProduct") << "Gen jets not found. Branch will not be filled" << std::endl;
134  }
135 
137  iEvent.getByToken(genParticleToken_, genParticles);
138 
139  if (genParticles.isValid()) {
140  int nPart{0};
141 
142  for (size_t i = 0; i < genParticles->size(); ++i) {
143  const reco::GenParticle& p = (*genParticles)[i];
144  int id = p.pdgId();
145 
146  // See if the parent was interesting
147  int parentID = -10000;
148  unsigned int nMo = p.numberOfMothers();
149  for (unsigned int i = 0; i < nMo; ++i) {
150  int thisParentID = dynamic_cast<const reco::GenParticle*>(p.mother(i))->pdgId();
151  //
152  // Is this a bottom hadron?
153  int hundredsIndex = abs(thisParentID) / 100;
154  int thousandsIndex = abs(thisParentID) / 1000;
155  if (((abs(thisParentID) >= 23) && (abs(thisParentID) <= 25)) || (abs(thisParentID) == 6) ||
156  (hundredsIndex == 5) || (hundredsIndex == 4) || (thousandsIndex == 5) || (thousandsIndex == 4))
157  parentID = thisParentID;
158  }
159  if ((parentID == -10000) && (nMo > 0))
160  parentID = dynamic_cast<const reco::GenParticle*>(p.mother(0))->pdgId();
161  //
162  // If the parent of this particle is interesting, store all of the info
163  if ((parentID != p.pdgId()) &&
164  ((parentID > -9999) || (abs(id) == 11) || (abs(id) == 13) || (abs(id) == 23) || (abs(id) == 24) ||
165  (abs(id) == 25) || (abs(id) == 4) || (abs(id) == 5) || (abs(id) == 6))) {
166  l1GenData_->partId.push_back(p.pdgId());
167  l1GenData_->partStat.push_back(p.status());
168  l1GenData_->partPt.push_back(p.pt());
169  l1GenData_->partEta.push_back(p.eta());
170  l1GenData_->partPhi.push_back(p.phi());
171  l1GenData_->partE.push_back(p.energy());
172  l1GenData_->partParent.push_back(parentID);
173  l1GenData_->partCh.push_back(p.charge());
174  ++nPart;
175  }
176  }
177  l1GenData_->nPart = nPart;
178  }
179 
181  iEvent.getByToken(pileupInfoToken_, puInfoCollection);
182 
183  if (!puInfoCollection.isValid()) {
184  throw cms::Exception("ProductNotValid") << "pileupInfoSource not valid";
185  }
186 
187  // Loop over vector, find in-time entry, then store the relevant info
188  std::vector<PileupSummaryInfo>::const_iterator puItr = puInfoCollection->begin();
189  std::vector<PileupSummaryInfo>::const_iterator puEnd = puInfoCollection->end();
190  for (; puItr != puEnd; ++puItr) {
191  int bx = puItr->getBunchCrossing();
192  if (bx == 0) {
193  l1GenData_->nMeanPU = puItr->getTrueNumInteractions();
194  l1GenData_->nVtx = puItr->getPU_NumInteractions();
195  break;
196  }
197  }
198 
199  tree_->Fill();
200 }
201 
202 // ------------ method called once each job just before starting event loop ------------
204 
205 // ------------ method called once each job just after ending the event loop ------------
207 
208 //define this as a plug-in
T getUntrackedParameter(std::string const &, T const &) const
const Candidate * mother(size_type=0) const override
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
double pt() const final
transverse momentum
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > pileupInfoToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void analyze(const edm::Event &, const edm::EventSetup &) override
int status() const final
status word
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::EDGetTokenT< GenEventInfoProduct > genInfoToken_
size_t numberOfMothers() const override
number of mothers
int pdgId() const final
PDG identifier.
~L1GenTreeProducer() override
edm::Service< TFileService > fs_
int iEvent
Definition: GenABIO.cc:224
void beginJob(void) override
edm::EDGetTokenT< reco::GenJetCollection > genJetToken_
edm::EDGetTokenT< reco::GenParticleCollection > genParticleToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:70
TString nPart(Int_t part, TString string, TString delimit=";", Bool_t removerest=true)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::unique_ptr< L1Analysis::L1AnalysisGeneratorDataFormat > l1GenData_
L1GenTreeProducer(const edm::ParameterSet &)
Log< level::Warning, false > LogWarning
void endJob() override
double phi() const final
momentum azimuthal angle
int charge() const final
electric charge
double energy() const final
energy
double eta() const final
momentum pseudorapidity