CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauCleaner.cc
Go to the documentation of this file.
1 /*
2  * RecoTauCleaner
3  *
4  * Author: Evan K. Friis, UC Davis
5  *
6  * Given a input collection of PFTaus, produces an output collection of PFTaus
7  * such that no two PFTaus come from the same PFJet. If multiple taus in the
8  * collection come from the same PFJet, (dirty) they are ordered according to a
9  * list of cleaners. Each cleaner is a RecoTauCleanerPlugin, and returns a
10  * double corresponding to the 'quality' for a given tau - an example would be
11  * the level of isolation. The set of dirty taus is then ranked
12  * lexicographically by these cleaners, and the best one is placed in the
13  * output collection.
14  *
15  * $Id $
16  */
17 
18 #include <boost/ptr_container/ptr_vector.hpp>
19 #include <algorithm>
20 #include <memory>
21 
27 
31 
34 
36 
37 template<typename Prod>
40  typedef boost::ptr_vector<Cleaner> CleanerList;
41  // Define our output type - i.e. reco::PFTau OR reco::PFTauRef
42  typedef typename Prod::value_type output_type;
43 
44  // Predicate that determines if two taus 'overlap' i.e. share a base PFJet
46  public:
48  { return a->jetRef() == b->jetRef(); }
49  };
50 
51  public:
52  explicit RecoTauCleanerImpl(const edm::ParameterSet& pset);
54  void produce(edm::Event& evt, const edm::EventSetup& es);
55 
56  private:
59  // Optional selection on the output of the taus
60  std::auto_ptr<StringCutObjectSelector<reco::PFTau> > outputSelector_;
61 };
62 
63 template<typename Prod>
65  tauSrc_ = pset.getParameter<edm::InputTag>("src");
66  // Build our list of quality plugins
67  typedef std::vector<edm::ParameterSet> VPSet;
68  // Get each of our tau builders
69  const VPSet& cleaners = pset.getParameter<VPSet>("cleaners");
70  for (VPSet::const_iterator cleanerPSet = cleaners.begin();
71  cleanerPSet != cleaners.end(); ++cleanerPSet) {
72  // Get plugin name
73  const std::string& pluginType =
74  cleanerPSet->getParameter<std::string>("plugin");
75  // Build the plugin
76  cleaners_.push_back(
77  RecoTauCleanerPluginFactory::get()->create(pluginType, *cleanerPSet));
78  }
79 
80  // Check if we want to apply a final output selection
81  if (pset.exists("outputSelection")) {
82  std::string selection = pset.getParameter<std::string>("outputSelection");
83  if (selection != "") {
84  outputSelector_.reset(
86  }
87  }
88 
89  // Build the predicate that ranks our taus.
90  produces<Prod>();
91 }
92 
93 namespace {
94 // Template to convert a ref to desired output type
95 template<typename T> const T convert(const reco::PFTauRef &tau);
96 
97 template<> const edm::RefToBase<reco::PFTau>
98 convert<edm::RefToBase<reco::PFTau> >(const reco::PFTauRef &tau) {
100 }
101 
102 template<> const reco::PFTauRef
103 convert<reco::PFTauRef>(const reco::PFTauRef &tau) { return tau; }
104 
105 template<> const reco::PFTau
106 convert<reco::PFTau>(const reco::PFTauRef &tau) { return *tau; }
107 }
108 
109 template<typename Prod>
111  const edm::EventSetup& es) {
112  // Update all our cleaners with the event info if they need it
113  for (CleanerList::iterator cleaner = cleaners_.begin();
114  cleaner != cleaners_.end(); ++cleaner) {
115  cleaner->setup(evt, es);
116  }
117 
118  // Get the input collection to clean
120  evt.getByLabel(tauSrc_, input);
121 
122  // Cast the input candidates to Refs to real taus
123  reco::PFTauRefVector inputRefs =
124  reco::tau::castView<reco::PFTauRefVector>(input);
125 
126  // Make an STL algorithm friendly vector of refs
127  typedef std::vector<reco::PFTauRef> PFTauRefs;
128  // Collection of all taus. Some are from the same PFJet. We must clean them.
129  PFTauRefs dirty(inputRefs.size());
130 
131  // Sort the input tau refs according to our predicate
132  auto N = inputRefs.size();
133  auto CN = cleaners_.size();
134  int index[N];
135  for (decltype(N) i=0; i!=N; ++i) index[i]=i;
136  float ranks[N*CN];
137  // auto ranks = [rr,CN](int i, int j){return rr[i*CN+j];};
138  // float const * * rr = ranks;
139  decltype(N) ii=0;
140  for (decltype(N) ir=0; ir!=N; ++ir)
141  for(decltype(N) cl=0;cl!=CN;++cl)
142  ranks[ii++]=cleaners_[cl](inputRefs[ir]);
143  assert(ii==(N*CN));
144  const float * rr = ranks;
145  std::sort(index,index+N,[rr,CN](int i, int j) { return std::lexicographical_compare(rr+i*CN,rr+i*CN+CN,rr+j*CN,rr+j*CN+CN);});
146  for (decltype(N) ir=0; ir!=N; ++ir)
147  dirty[ir]=inputRefs[index[ir]];
148 
149  // Clean the taus, ensuring that only one tau per jet is produced
150  PFTauRefs cleanTaus = reco::tau::cleanOverlaps<PFTauRefs,
151  RemoveDuplicateJets>(dirty);
152 
153  // create output collection
154  std::auto_ptr<Prod> output(new Prod());
155  //output->reserve(cleanTaus.size());
156 
157  // Copy clean refs into output
158  for (PFTauRefs::const_iterator tau = cleanTaus.begin();
159  tau != cleanTaus.end(); ++tau) {
160  // If we are applying an output selection, check if it passes
161  bool selected = true;
162  if (outputSelector_.get() && !(*outputSelector_)(**tau)) {
163  selected = false;
164  }
165  if (selected) {
166  output->push_back(convert<output_type>(*tau));
167  }
168  }
169  evt.put(output);
170 }
171 
174 
178 
void produce(edm::Event &evt, const edm::EventSetup &es)
T getParameter(std::string const &) const
RecoTauCleanerImpl< reco::PFTauRefVector > RecoTauRefCleaner
int i
Definition: DBlmapReader.cc:9
static float CN[]
Definition: sicif.h:76
edm::InputTag tauSrc_
std::auto_ptr< StringCutObjectSelector< reco::PFTau > > outputSelector_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool operator()(const reco::PFTauRef &a, const reco::PFTauRef &b)
selection
main part
Definition: corrVsCorr.py:98
bool exists(std::string const &parameterName) const
checks if a parameter exists
RecoTauCleanerImpl< reco::PFTauCollection > RecoTauCleaner
int ii
Definition: cuy.py:588
void convert(uint32 i, char_uint32 v)
Definition: MsgTools.h:46
CleanerList cleaners_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
int j
Definition: DBlmapReader.cc:9
Container::value_type value_type
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
#define N
Definition: blowfish.cc:9
double b
Definition: hdecay.h:120
reco::tau::RecoTauCleanerPlugin Cleaner
double a
Definition: hdecay.h:121
boost::ptr_vector< Cleaner > CleanerList
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
Container cleanOverlaps(const Container &dirty)
long double T
Prod::value_type output_type
SurfaceDeformation * create(int type, const std::vector< double > &params)
T get(const Candidate &c)
Definition: component.h:56
RecoTauCleanerImpl(const edm::ParameterSet &pset)