CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauProducer.cc
Go to the documentation of this file.
1 /*
2  * RecoTauProducer
3  *
4  * Interface between the various tau algorithms and the edm::Event. The
5  * RecoTauProducer takes as data input is a collection (view) of reco::PFJets,
6  * and Jet-PiZero assoications that give the reco::RecoTauPiZeros for those
7  * jets. The actaul building of taus is done by the list of builders - each of
8  * which constructs a PFTau for each PFJet. The output collection may have
9  * multiple taus for each PFJet - these overlaps are to be resolved by the
10  * RecoTauCleaner module.
11  *
12  * Additionally, there are "modifier" plugins, which can do things like add the
13  * lead track significance, or electron rejection variables.
14  *
15  * Author: Evan K. Friis (UC Davis)
16  *
17  */
18 #include <boost/ptr_container/ptr_vector.hpp>
19 #include <boost/foreach.hpp>
20 
21 #include <algorithm>
22 #include <functional>
23 
29 
32 
37 
39 
41  public:
44  typedef boost::ptr_vector<Builder> BuilderList;
45  typedef boost::ptr_vector<Modifier> ModifierList;
46 
47  explicit RecoTauProducer(const edm::ParameterSet& pset);
49  void produce(edm::Event& evt, const edm::EventSetup& es);
50 
51  private:
56  // Optional selection on the output of the taus
57  std::auto_ptr<StringCutObjectSelector<reco::PFTau> > outputSelector_;
58  // Whether or not to add build a tau from a jet for which the builders
59  // return no taus. The tau will have no content, only the four vector of
60  // the orginal jet.
62 };
63 
65  jetSrc_ = pset.getParameter<edm::InputTag>("jetSrc");
66  piZeroSrc_ = pset.getParameter<edm::InputTag>("piZeroSrc");
67 
68  typedef std::vector<edm::ParameterSet> VPSet;
69  // Get each of our tau builders
70  const VPSet& builders = pset.getParameter<VPSet>("builders");
71  for (VPSet::const_iterator builderPSet = builders.begin();
72  builderPSet != builders.end(); ++builderPSet) {
73  // Get plugin name
74  const std::string& pluginType =
75  builderPSet->getParameter<std::string>("plugin");
76  // Build the plugin
77  builders_.push_back(
78  RecoTauBuilderPluginFactory::get()->create(pluginType, *builderPSet));
79  }
80 
81  const VPSet& modfiers = pset.getParameter<VPSet>("modifiers");
82  for (VPSet::const_iterator modfierPSet = modfiers.begin();
83  modfierPSet != modfiers.end(); ++modfierPSet) {
84  // Get plugin name
85  const std::string& pluginType =
86  modfierPSet->getParameter<std::string>("plugin");
87  // Build the plugin
88  modifiers_.push_back(
89  RecoTauModifierPluginFactory::get()->create(pluginType, *modfierPSet));
90  }
91 
92  // Check if we want to apply a final output selection
93  if (pset.exists("outputSelection")) {
94  std::string selection = pset.getParameter<std::string>("outputSelection");
95  if (selection != "") {
96  outputSelector_.reset(
98  }
99  }
100  buildNullTaus_ = pset.getParameter<bool>("buildNullTaus");
101  produces<reco::PFTauCollection>();
102 }
103 
105  // Get the jet input collection via a view of Candidates
107  evt.getByLabel(jetSrc_, jetView);
108 
109  // Convert to a vector of PFJetRefs
111  reco::tau::castView<reco::PFJetRefVector>(jetView);
112 
113  // Get the pizero input collection
115  evt.getByLabel(piZeroSrc_, piZeroAssoc);
116 
117  // Update all our builders and modifiers with the event info
118  for (BuilderList::iterator builder = builders_.begin();
119  builder != builders_.end(); ++builder) {
120  builder->setup(evt, es);
121  }
122  for (ModifierList::iterator modifier = modifiers_.begin();
123  modifier != modifiers_.end(); ++modifier) {
124  modifier->setup(evt, es);
125  }
126 
127  // Create output collection
128  std::auto_ptr<reco::PFTauCollection> output(new reco::PFTauCollection());
129 
130  // Loop over the jets and build the taus for each jet
131  BOOST_FOREACH(reco::PFJetRef jetRef, jets) {
132  // Get the PiZeros associated with this jet
133  const std::vector<reco::RecoTauPiZero>& piZeros = (*piZeroAssoc)[jetRef];
134  // Loop over our builders and create the set of taus for this jet
135  unsigned int nTausBuilt = 0;
136  for (BuilderList::const_iterator builder = builders_.begin();
137  builder != builders_.end(); ++builder) {
138  // Get a ptr_vector of taus from the builder
140  (*builder)(jetRef, piZeros));
141  // Check this size of the taus built.
142  // Grow the vector if necessary
143  output->reserve(output->size() + taus.size());
144  // Copy without selection
145  if (!outputSelector_.get()) {
146  output->insert(output->end(), taus.begin(), taus.end());
147  nTausBuilt += taus.size();
148  } else {
149  // Copy only those that pass the selection.
150  BOOST_FOREACH(const reco::PFTau& tau, taus) {
151  if ((*outputSelector_)(tau)) {
152  nTausBuilt++;
153  output->push_back(tau);
154  }
155  }
156  }
157  }
158  // If we didn't build *any* taus for this jet, build a null tau if desired.
159  // The null PFTau has no content, but it's four vector is set to that of the
160  // jet.
161  if (!nTausBuilt && buildNullTaus_) {
162  reco::PFTau nullTau(std::numeric_limits<int>::quiet_NaN(), jetRef->p4());
163  nullTau.setjetRef(jetRef);
164  output->push_back(nullTau);
165  }
166  }
167 
168  // Loop over the taus we have created and apply our modifiers to the taus
169  for (reco::PFTauCollection::iterator tau = output->begin();
170  tau != output->end(); ++tau) {
171  for (ModifierList::const_iterator modifier = modifiers_.begin();
172  modifier != modifiers_.end(); ++modifier) {
173  (*modifier)(*tau);
174  }
175  }
176  evt.put(output);
177 }
178 
T getParameter(std::string const &) const
edm::InputTag piZeroSrc_
boost::ptr_vector< Modifier > ModifierList
std::vector< PFTau > PFTauCollection
collection of PFTau objects
Definition: PFTauFwd.h:9
BuilderList builders_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
boost::ptr_vector< reco::PFTau > output_type
ModifierList modifiers_
std::auto_ptr< StringCutObjectSelector< reco::PFTau > > outputSelector_
void setjetRef(const PFJetRef &)
Definition: PFTau.cc:51
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
RecoTauProducer(const edm::ParameterSet &pset)
tuple pset
Definition: CrabTask.py:85
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
edm::InputTag jetSrc_
reco::tau::RecoTauModifierPlugin Modifier
reco::tau::RecoTauBuilderPlugin Builder
void produce(edm::Event &evt, const edm::EventSetup &es)
boost::ptr_vector< Builder > BuilderList
SurfaceDeformation * create(int type, const std::vector< double > &params)
T get(const Candidate &c)
Definition: component.h:56