CMS 3D CMS Logo

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 actual 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  * Authors: Evan K. Friis (UC Davis),
16  * Christian Veelken (LLR)
17  *
18  */
19 #include "boost/bind.hpp"
20 #include <boost/ptr_container/ptr_vector.hpp>
21 
22 #include <algorithm>
23 #include <functional>
24 
30 
33 
41 
43 
45 {
46  public:
49  typedef boost::ptr_vector<Builder> BuilderList;
50  typedef boost::ptr_vector<Modifier> ModifierList;
51 
52  explicit RecoTauProducer(const edm::ParameterSet& pset);
53  ~RecoTauProducer() override {}
54  void produce(edm::Event& evt, const edm::EventSetup& es) override;
55 
56  private:
61 
62  double minJetPt_;
63  double maxJetAbsEta_;
64  //token definition
69 
70  BuilderList builders_;
71  ModifierList modifiers_;
72  // Optional selection on the output of the taus
73  std::unique_ptr<StringCutObjectSelector<reco::PFTau> > outputSelector_;
74  // Whether or not to add build a tau from a jet for which the builders
75  // return no taus. The tau will have no content, only the four vector of
76  // the orginal jet.
78 };
79 
81 {
82  jetSrc_ = pset.getParameter<edm::InputTag>("jetSrc");
83  jetRegionSrc_ = pset.getParameter<edm::InputTag>("jetRegionSrc");
84  chargedHadronSrc_ = pset.getParameter<edm::InputTag>("chargedHadronSrc");
85  piZeroSrc_ = pset.getParameter<edm::InputTag>("piZeroSrc");
86 
87  minJetPt_ = ( pset.exists("minJetPt") ) ? pset.getParameter<double>("minJetPt") : -1.0;
88  maxJetAbsEta_ = ( pset.exists("maxJetAbsEta") ) ? pset.getParameter<double>("maxJetAbsEta") : 99.0;
89  //consumes definition
90  jet_token=consumes<reco::CandidateView>(jetSrc_);
91  jetRegion_token = consumes<edm::Association<reco::PFJetCollection> >(jetRegionSrc_);
92  chargedHadron_token = consumes<reco::PFJetChargedHadronAssociation>(chargedHadronSrc_);
93  piZero_token = consumes<reco::JetPiZeroAssociation>(piZeroSrc_);
94 
95  typedef std::vector<edm::ParameterSet> VPSet;
96  // Get each of our tau builders
97  const VPSet& builders = pset.getParameter<VPSet>("builders");
98  for ( VPSet::const_iterator builderPSet = builders.begin();
99  builderPSet != builders.end(); ++builderPSet ) {
100  // Get plugin name
101  const std::string& pluginType = builderPSet->getParameter<std::string>("plugin");
102  // Build the plugin
103  builders_.push_back(RecoTauBuilderPluginFactory::get()->create(pluginType, *builderPSet, consumesCollector()));
104  }
105 
106  const VPSet& modfiers = pset.getParameter<VPSet>("modifiers");
107  for ( VPSet::const_iterator modfierPSet = modfiers.begin();
108  modfierPSet != modfiers.end(); ++modfierPSet) {
109  // Get plugin name
110  const std::string& pluginType = modfierPSet->getParameter<std::string>("plugin");
111  // Build the plugin
113  plugin = RecoTauModifierPluginFactory::get()->create(pluginType, *modfierPSet, consumesCollector());
114  modifiers_.push_back(plugin);
115  }
116 
117  // Check if we want to apply a final output selection
118  if ( pset.exists("outputSelection") ) {
119  std::string selection = pset.getParameter<std::string>("outputSelection");
120  if ( !selection.empty() ) {
122  }
123  }
124  buildNullTaus_ = pset.getParameter<bool>("buildNullTaus");
125 
126  produces<reco::PFTauCollection>();
127 }
128 
130 {
131  // Get the jet input collection via a view of Candidates
133  evt.getByToken(jet_token, jetView);
134 
135  // Convert to a vector of PFJetRefs
136  reco::PFJetRefVector jets = reco::tau::castView<reco::PFJetRefVector>(jetView);
137 
138  // Get the jet region producer
140  evt.getByToken(jetRegion_token, jetRegionHandle);
141 
142  // Get the charged hadron input collection
144  evt.getByToken(chargedHadron_token, chargedHadronAssoc);
145 
146  // Get the pizero input collection
148  evt.getByToken(piZero_token, piZeroAssoc);
149 
150  // Update all our builders and modifiers with the event info
151  for (BuilderList::iterator builder = builders_.begin();
152  builder != builders_.end(); ++builder) {
153  builder->setup(evt, es);
154  }
155  for (ModifierList::iterator modifier = modifiers_.begin();
156  modifier != modifiers_.end(); ++modifier) {
157  modifier->setup(evt, es);
158  }
159 
160  // Create output collection
161  auto output = std::make_unique<reco::PFTauCollection>();
162  output->reserve(jets.size());
163 
164  // Loop over the jets and build the taus for each jet
165  for(auto const& jetRef : jets ) {
166  // Get the jet with extra constituents from an area around the jet
167  if(jetRef->pt() - minJetPt_ < 1e-5) continue;
168  if(std::abs(jetRef->eta()) - maxJetAbsEta_ > -1e-5) continue;
169  reco::PFJetRef jetRegionRef = (*jetRegionHandle)[jetRef];
170  if ( jetRegionRef.isNull() ) {
171  throw cms::Exception("BadJetRegionRef")
172  << "No jet region can be found for the current jet: " << jetRef.id();
173  }
174  // Remove all the jet constituents from the jet extras
175  std::vector<reco::PFCandidatePtr> jetCands = jetRef->getPFConstituents();
176  std::vector<reco::PFCandidatePtr> allRegionalCands = jetRegionRef->getPFConstituents();
177  // Sort both by ref key
178  std::sort(jetCands.begin(), jetCands.end());
179  std::sort(allRegionalCands.begin(), allRegionalCands.end());
180  // Get the regional junk candidates not in the jet.
181  std::vector<reco::PFCandidatePtr> uniqueRegionalCands;
182 
183  // This can actually be less than zero, if the jet has really crazy soft
184  // stuff really far away from the jet axis.
185  if ( allRegionalCands.size() > jetCands.size() ) {
186  uniqueRegionalCands.reserve(allRegionalCands.size() - jetCands.size());
187  }
188 
189  // Subtract the jet cands from the regional cands
190  std::set_difference(allRegionalCands.begin(), allRegionalCands.end(),
191  jetCands.begin(), jetCands.end(),
192  std::back_inserter(uniqueRegionalCands));
193 
194  // Get the charged hadrons associated with this jet
195  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons = (*chargedHadronAssoc)[jetRef];
196 
197  // Get the pizeros associated with this jet
198  const std::vector<reco::RecoTauPiZero>& piZeros = (*piZeroAssoc)[jetRef];
199  // Loop over our builders and create the set of taus for this jet
200  unsigned int nTausBuilt = 0;
201  for ( BuilderList::const_iterator builder = builders_.begin();
202  builder != builders_.end(); ++builder) {
203  // Get a ptr_vector of taus from the builder
204  reco::tau::RecoTauBuilderPlugin::output_type taus((*builder)(jetRef, chargedHadrons, piZeros, uniqueRegionalCands));
205  // Make sure all taus have their jetref set correctly
206  std::for_each(taus.begin(), taus.end(), boost::bind(&reco::PFTau::setjetRef, _1, jetRef));
207  // Copy without selection
208  if ( !outputSelector_.get() ) {
209  output->insert(output->end(), taus.begin(), taus.end());
210  nTausBuilt += taus.size();
211  } else {
212  // Copy only those that pass the selection.
213  for(auto const& tau : taus ) {
214  if ( (*outputSelector_)(tau) ) {
215  nTausBuilt++;
216  output->push_back(tau);
217  }
218  }
219  }
220  }
221  // If we didn't build *any* taus for this jet, build a null tau if desired.
222  // The null PFTau has no content, but it's four vector is set to that of the
223  // jet.
224  if ( !nTausBuilt && buildNullTaus_ ) {
225  reco::PFTau nullTau(std::numeric_limits<int>::quiet_NaN(), jetRef->p4());
226  nullTau.setjetRef(jetRef);
227  output->push_back(nullTau);
228  }
229  }
230 
231  // Loop over the taus we have created and apply our modifiers to the taus
232  for ( reco::PFTauCollection::iterator tau = output->begin();
233  tau != output->end(); ++tau ) {
234  for ( ModifierList::const_iterator modifier = modifiers_.begin();
235  modifier != modifiers_.end(); ++modifier ) {
236  (*modifier)(*tau);
237  }
238  }
239 
240  for ( ModifierList::iterator modifier = modifiers_.begin();
241  modifier != modifiers_.end(); ++modifier ) {
242  modifier->endEvent();
243  }
244 
245  evt.put(std::move(output));
246 }
247 
T getParameter(std::string const &) const
edm::InputTag piZeroSrc_
boost::ptr_vector< Modifier > ModifierList
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
edm::EDGetTokenT< reco::CandidateView > jet_token
~RecoTauProducer() override
def create(alignables, pedeDump, additionalData, outputFile, config)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
BuilderList builders_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
selection
main part
Definition: corrVsCorr.py:99
bool exists(std::string const &parameterName) const
checks if a parameter exists
boost::ptr_vector< reco::PFTau > output_type
ModifierList modifiers_
void setjetRef(const PFJetRef &)
Definition: PFTau.cc:59
RecoTauProducer(const edm::ParameterSet &pset)
unique_ptr< JetDefinition::Plugin > plugin
vector< PseudoJet > jets
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::InputTag chargedHadronSrc_
bool isNull() const
Checks for null.
Definition: Ref.h:248
edm::InputTag jetSrc_
void produce(edm::Event &evt, const edm::EventSetup &es) override
edm::EDGetTokenT< reco::JetPiZeroAssociation > piZero_token
reco::tau::RecoTauModifierPlugin Modifier
edm::EDGetTokenT< reco::PFJetChargedHadronAssociation > chargedHadron_token
edm::EDGetTokenT< edm::Association< reco::PFJetCollection > > jetRegion_token
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
reco::tau::RecoTauBuilderPlugin Builder
edm::InputTag jetRegionSrc_
boost::ptr_vector< Builder > BuilderList
def move(src, dest)
Definition: eostools.py:511
T get(const Candidate &c)
Definition: component.h:55
std::unique_ptr< StringCutObjectSelector< reco::PFTau > > outputSelector_