CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TagProbeFitTreeProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: TagProbeFitTreeProducer
4 // Class: TagProbeFitTreeProducer
5 //
13 //
14 // Original Author: Sep 15 09:45
15 // Created: Mon Sep 15 09:49:08 CEST 2008
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 #include <ctype.h>
23 
24 // user include files
27 
32 
35 
38 
41 
42 #include <set>
43 
44 //
45 // class decleration
46 //
47 
49  public:
52 
53 
54  private:
55  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
56  virtual void endJob() override ;
57 
58  //---- MC truth information
60  bool isMC_;
64  std::set<int32_t> motherPdgId_;
66  bool checkMother(const reco::GenParticleRef &ref) const ;
67 
68  //---- Unbiased MC truth information
75 
79  std::auto_ptr<tnp::TPTreeFiller> treeFiller_;
81  std::auto_ptr<tnp::BaseTreeFiller> mcUnbiasFiller_;
82  std::auto_ptr<tnp::BaseTreeFiller> oldTagFiller_;
83  std::auto_ptr<tnp::BaseTreeFiller> tagFiller_;
84  std::auto_ptr<tnp::BaseTreeFiller> pairFiller_;
85  std::auto_ptr<tnp::BaseTreeFiller> mcFiller_;
86 };
87 
88 //
89 // constructors and destructor
90 //
92  isMC_(iConfig.getParameter<bool>("isMC")),
93  makeMCUnbiasTree_(isMC_ ? iConfig.getParameter<bool>("makeMCUnbiasTree") : false),
94  checkMotherInUnbiasEff_(makeMCUnbiasTree_ ? iConfig.getParameter<bool>("checkMotherInUnbiasEff") : false),
95  tagProbePairMaker_(iConfig, consumesCollector()),
96  treeFiller_(new tnp::TPTreeFiller(iConfig, consumesCollector())),
97  oldTagFiller_((iConfig.existsAs<bool>("fillTagTree") && iConfig.getParameter<bool>("fillTagTree")) ? new tnp::BaseTreeFiller("tag_tree",iConfig, consumesCollector()) : 0)
98 {
99  if (isMC_) {
100  // For mc efficiency we need the MC matches for tags & probes
101  tagMatchesToken_ = consumes<edm::Association<std::vector<reco::GenParticle> > >(iConfig.getParameter<edm::InputTag>("tagMatches"));
102  probeMatchesToken_ = consumes<edm::Association<std::vector<reco::GenParticle> > >(iConfig.getParameter<edm::InputTag>("probeMatches"));
103  //.. and the pdgids of the possible mothers
104  if (iConfig.existsAs<int32_t>("motherPdgId")) {
105  motherPdgId_.insert(iConfig.getParameter<int32_t>("motherPdgId"));
106  } else {
107  std::vector<int32_t> motherIds = iConfig.getParameter<std::vector<int32_t> >("motherPdgId");
108  motherPdgId_.insert(motherIds.begin(), motherIds.end());
109  }
110 
111  // For unbiased efficiency we also need the collection of all probes
112  if (makeMCUnbiasTree_) {
113  allProbesToken_ = consumes<reco::CandidateView>(iConfig.getParameter<edm::InputTag>("allProbes"));
114  mcUnbiasFiller_.reset(new tnp::BaseTreeFiller("mcUnbias_tree",iConfig, consumesCollector()));
115  }
116  }
117 
118  edm::ParameterSet tagPSet;
119  if (iConfig.existsAs<edm::ParameterSet>("tagVariables")) tagPSet.addParameter<edm::ParameterSet>("variables", iConfig.getParameter<edm::ParameterSet>("tagVariables"));
120  if (iConfig.existsAs<edm::ParameterSet>("tagFlags" )) tagPSet.addParameter<edm::ParameterSet>("flags", iConfig.getParameter<edm::ParameterSet>("tagFlags"));
121  if (!tagPSet.empty()) { tagFiller_.reset(new tnp::BaseTreeFiller(*treeFiller_, tagPSet, consumesCollector(), "tag_")); }
122  edm::ParameterSet mcPSet;
123  if (iConfig.existsAs<edm::ParameterSet>("mcVariables")) mcPSet.addParameter<edm::ParameterSet>("variables", iConfig.getParameter<edm::ParameterSet>("mcVariables"));
124  if (iConfig.existsAs<edm::ParameterSet>("mcFlags" )) mcPSet.addParameter<edm::ParameterSet>("flags", iConfig.getParameter<edm::ParameterSet>("mcFlags"));
125  if (!mcPSet.empty()) { mcFiller_.reset(new tnp::BaseTreeFiller(*treeFiller_, mcPSet, consumesCollector(), "mc_")); }
126  edm::ParameterSet pairPSet;
127  if (iConfig.existsAs<edm::ParameterSet>("pairVariables")) pairPSet.addParameter<edm::ParameterSet>("variables", iConfig.getParameter<edm::ParameterSet>("pairVariables"));
128  if (iConfig.existsAs<edm::ParameterSet>("pairFlags" )) pairPSet.addParameter<edm::ParameterSet>("flags", iConfig.getParameter<edm::ParameterSet>("pairFlags"));
129  if (!pairPSet.empty()) { pairFiller_.reset(new tnp::BaseTreeFiller(*treeFiller_, pairPSet, consumesCollector(), "pair_")); }
130 }
131 
132 
134 {
135 }
136 
137 
138 //
139 // member functions
140 //
141 
142 // ------------ method called to for each event ------------
143 void
145 {
146  using namespace edm; using namespace std;
147  Handle<reco::CandidateView> src, allProbes;
148  Handle<Association<vector<reco::GenParticle> > > tagMatches, probeMatches;
149  treeFiller_->init(iEvent); // read out info from the event if needed (external vars, list of passing probes, ...)
150  if (oldTagFiller_.get()) oldTagFiller_->init(iEvent);
151  if (tagFiller_.get()) tagFiller_->init(iEvent);
152  if (pairFiller_.get()) pairFiller_->init(iEvent);
153  if (mcFiller_.get()) mcFiller_->init(iEvent);
154 
155  // on mc we want to load also the MC match info
156  if (isMC_) {
157  iEvent.getByToken(tagMatchesToken_, tagMatches);
158  iEvent.getByToken(probeMatchesToken_, probeMatches);
159  }
160 
161  // get the list of (tag+probe) pairs, performing arbitration
162  tnp::TagProbePairs pairs = tagProbePairMaker_.run(iEvent);
163  // loop on them to fill the tree
164  for (tnp::TagProbePairs::const_iterator it = pairs.begin(), ed = pairs.end(); it != ed; ++it) {
165  // on mc, fill mc info (on non-mc, let it to 'true', the treeFiller will ignore it anyway
166  bool mcTrue = false;
167  float mcMass = 0.f;
168  if (isMC_) {
169  reco::GenParticleRef mtag = (*tagMatches)[it->tag], mprobe = (*probeMatches)[it->probe];
170  mcTrue = checkMother(mtag) && checkMother(mprobe);
171  if (mcTrue) {
172  mcMass = (mtag->p4() + mprobe->p4()).mass();
173  if (mcFiller_.get()) mcFiller_->fill(reco::CandidateBaseRef(mprobe));
174  }
175  }
176  // fill in the variables for this t+p pair
177  if (tagFiller_.get()) tagFiller_->fill(it->tag);
178  if (oldTagFiller_.get()) oldTagFiller_->fill(it->tag);
179  if (pairFiller_.get()) pairFiller_->fill(it->pair);
180  treeFiller_->fill(it->probe, it->mass, mcTrue, mcMass);
181  }
182 
183  if (isMC_ && makeMCUnbiasTree_) {
184  // read full collection of probes
185  iEvent.getByToken(allProbesToken_, allProbes);
186  // init the tree filler
187  mcUnbiasFiller_->init(iEvent);
188  // loop on probes
189  for (size_t i = 0, n = allProbes->size(); i < n; ++i) {
190  const reco::CandidateBaseRef & probe = allProbes->refAt(i);
191  // check mc match, and possibly mother match
192  reco::GenParticleRef probeMatch = (*probeMatches)[probe];
193  bool probeOk = checkMotherInUnbiasEff_ ? checkMother(probeMatch) : probeMatch.isNonnull();
194  // fill the tree only for good ones
195  if (probeOk) mcUnbiasFiller_->fill(probe);
196  }
197  }
198 
199 }
200 
201 bool
203  if (ref.isNull()) return false;
204  if (motherPdgId_.empty()) return true;
205  if (motherPdgId_.find(abs(ref->pdgId())) != motherPdgId_.end()) return true;
206  reco::GenParticle::mothers m = ref->motherRefVector();
207  for (reco::GenParticle::mothers::const_iterator it = m.begin(), e = m.end(); it != e; ++it) {
208  if (checkMother(*it)) return true;
209  }
210  return false;
211 }
212 
213 
214 // ------------ method called once each job just after ending the event loop ------------
215 void
217  // ask to write the current PSet info into the TTree header
219 }
220 
221 //define this as a plug-in
T getParameter(std::string const &) const
bool empty() const
Definition: ParameterSet.h:218
edm::EDGetTokenT< edm::Association< std::vector< reco::GenParticle > > > tagMatchesToken_
Token foran edm::Association&lt;reco::GenParticle&gt; from tags &amp; probes to MC truth.
int i
Definition: DBlmapReader.cc:9
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
std::auto_ptr< tnp::TPTreeFiller > treeFiller_
The object that actually computes variables and fills the tree for T&amp;P.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
std::vector< TagProbePair > TagProbePairs
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::auto_ptr< tnp::BaseTreeFiller > oldTagFiller_
std::auto_ptr< tnp::BaseTreeFiller > tagFiller_
edm::EDGetTokenT< reco::CandidateView > allProbesToken_
InputTag to the collection of all probes.
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:255
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:250
ParameterSet const & getProcessParameterSetContainingModule(ModuleDescription const &moduleDescription)
bool checkMotherInUnbiasEff_
Check mother pdgId in unbiased inefficiency measurement.
ModuleDescription const & moduleDescription() const
Definition: EDAnalyzer.h:43
int iEvent
Definition: GenABIO.cc:230
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
bool isMC_
Is this sample MC?
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:144
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::auto_ptr< tnp::BaseTreeFiller > mcFiller_
virtual void endJob() override
bool makeMCUnbiasTree_
Do we have to compute this.
bool isNull() const
Checks for null.
Definition: Ref.h:249
tnp::TagProbePairMaker tagProbePairMaker_
The object that produces pairs of tags and probes, making any arbitration needed. ...
std::auto_ptr< tnp::BaseTreeFiller > mcUnbiasFiller_
The object that actually computes variables and fills the tree for unbiased MC.
edm::EDGetTokenT< edm::Association< std::vector< reco::GenParticle > > > probeMatchesToken_
bool checkMother(const reco::GenParticleRef &ref) const
Return true if ref is not null and has an ancestor with pdgId inside &#39;motherPdgId_&#39;.
std::set< int32_t > motherPdgId_
Possible pdgids for the mother. If empty, any truth-matched mu will be considered good...
TagProbePairs run(const edm::Event &iEvent) const
fill in tghe T&amp;P pairs for this event
volatile std::atomic< bool > shutdown_flag false
std::auto_ptr< tnp::BaseTreeFiller > pairFiller_
TagProbeFitTreeProducer(const edm::ParameterSet &)