CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauBuilderCombinatoricPlugin.cc
Go to the documentation of this file.
1 #include <vector>
2 
5 
12 
13 namespace reco { namespace tau {
14 
16  public:
20  const std::vector<RecoTauPiZero>& piZeros) const;
21  private:
25  struct decayModeInfo {
26  uint32_t maxPiZeros_;
27  uint32_t maxPFCHs_;
28  uint32_t nCharged_;
29  uint32_t nPiZeros_;
30  };
31  std::vector<decayModeInfo> decayModesToBuild_;
32 };
33 
36  qcuts_(pset.getParameter<edm::ParameterSet>("qualityCuts")),
37  usePFLeptonsAsChargedHadrons_(pset.getParameter<bool>("usePFLeptons")),
38  isolationConeSize_(pset.getParameter<double>("isolationConeSize")) {
39  typedef std::vector<edm::ParameterSet> VPSet;
40  const VPSet& decayModes = pset.getParameter<VPSet>("decayModes");
41  for (VPSet::const_iterator dm = decayModes.begin();
42  dm != decayModes.end(); ++dm) {
44  info.nCharged_ = dm->getParameter<uint32_t>("nCharged");
45  info.nPiZeros_ = dm->getParameter<uint32_t>("nPiZeros");
46  info.maxPFCHs_ = dm->getParameter<uint32_t>("maxTracks");
47  info.maxPiZeros_ = dm->getParameter<uint32_t>("maxPiZeros");
48  decayModesToBuild_.push_back(info);
49  }
50 }
51 
54  const reco::PFJetRef& jet,
55  const std::vector<RecoTauPiZero>& piZeros) const {
56  typedef std::vector<PFCandidatePtr> PFCandPtrs;
57  typedef std::vector<RecoTauPiZero> PiZeroList;
58 
60 
61  // Update the primary vertex used by the quality cuts. The PV is supplied by
62  // the base class.
64 
65  // Get PFCHs from this jet. They are already sorted by descending Pt
66  PFCandPtrs pfchs;
69  } else {
70  // Check if we want to include electrons in muons in "charged hadron"
71  // collection. This is the preferred behavior, as the PF lepton selections
72  // are very loose.
73  pfchs = qcuts_.filterRefs(pfChargedCands(*jet));
74  }
75 
76  PFCandPtrs pfnhs = qcuts_.filterRefs(
78 
79  // Loop over the decay modes we want to build
80  for (std::vector<decayModeInfo>::const_iterator
81  decayMode = decayModesToBuild_.begin();
82  decayMode != decayModesToBuild_.end(); ++decayMode) {
83  // Find how many piZeros are in this decay mode
84  size_t piZerosToBuild = decayMode->nPiZeros_;
85  // Find how many tracks are in this decay mode
86  size_t tracksToBuild = decayMode->nCharged_;
87 
88  // Skip decay mode if jet doesn't have the multiplicity to support it
89  if (pfchs.size() < tracksToBuild || piZeros.size() < piZerosToBuild)
90  continue;
91 
92  // Find the start and end of potential signal tracks
93  PFCandPtrs::iterator pfch_begin = pfchs.begin();
94  PFCandPtrs::iterator pfch_end = pfchs.end();
95  pfch_end = takeNElements(pfch_begin, pfch_end, decayMode->maxPFCHs_);
96 
97  // Build our track combo generator
99  PFCombo trackCombos(pfch_begin, pfch_end, tracksToBuild);
100 
101  // Find the start and end of potential signal tracks
102  PiZeroList::const_iterator piZero_begin = piZeros.begin();
103  PiZeroList::const_iterator piZero_end = piZeros.end();
104  piZero_end = takeNElements(piZero_begin, piZero_end,
105  decayMode->maxPiZeros_);
106 
107  // Build our piZero combo generator
108  typedef tau::CombinatoricGenerator<PiZeroList> PiZeroCombo;
109  PiZeroCombo piZeroCombos(piZero_begin, piZero_end, piZerosToBuild);
110 
111  /*
112  * Begin combinatoric loop for this decay mode
113  */
114 
115  // Loop over the different combinations of tracks
116  for (PFCombo::iterator trackCombo = trackCombos.begin();
117  trackCombo != trackCombos.end(); ++trackCombo) {
118  // Loop over the different combinations of PiZeros
119  for (PiZeroCombo::iterator piZeroCombo = piZeroCombos.begin();
120  piZeroCombo != piZeroCombos.end(); ++piZeroCombo) {
121  // Output tau
122  RecoTauConstructor tau(jet, getPFCands(), true);
123  // Reserve space in our collections
125  RecoTauConstructor::kChargedHadron, tracksToBuild);
126  tau.reserve(
128  RecoTauConstructor::kGamma, 2*piZerosToBuild); // k-factor = 2
129  tau.reservePiZero(RecoTauConstructor::kSignal, piZerosToBuild);
130 
131  // FIXME - are all these reserves okay? will they get propagated to the
132  // dataformat size if they are wrong?
133  tau.reserve(
135  RecoTauConstructor::kChargedHadron, pfchs.size() - tracksToBuild);
138  (piZeros.size() - piZerosToBuild)*2);
140  (piZeros.size() - piZerosToBuild));
141 
142  // Set signal and isolation components for charged hadrons, after
143  // converting them to a PFCandidateRefVector
144  tau.addPFCands(
146  trackCombo->combo_begin(), trackCombo->combo_end()
147  );
148 
149  // Get signal PiZero constituents and add them to the tau.
150  // The sub-gammas are automatically added.
151  tau.addPiZeros(
153  piZeroCombo->combo_begin(), piZeroCombo->combo_end()
154  );
155 
156  // Now build isolation collections
157  // Load our isolation tools
158  using namespace reco::tau::cone;
159  PFCandPtrDRFilter isolationConeFilter(tau.p4(), 0, isolationConeSize_);
160 
161  PiZeroDRFilter isolationConeFilterPiZero(
162  tau.p4(), 0, isolationConeSize_);
163 
164 
165  // Filter the isolation candidates in a DR cone
166  tau.addPFCands(
168  boost::make_filter_iterator(
169  isolationConeFilter,
170  trackCombo->remainder_begin(), trackCombo->remainder_end()),
171  boost::make_filter_iterator(
172  isolationConeFilter,
173  trackCombo->remainder_end(), trackCombo->remainder_end())
174  );
175 
176  // Add all the candidates that weren't included in the combinatoric
177  // generation
178  tau.addPFCands(
180  boost::make_filter_iterator(
181  isolationConeFilter,
182  pfch_end, pfchs.end()),
183  boost::make_filter_iterator(
184  isolationConeFilter,
185  pfchs.end(), pfchs.end())
186  );
187 
188  // Add all the netural hadron candidates to the isolation collection
189  tau.addPFCands(
191  boost::make_filter_iterator(
192  isolationConeFilter,
193  pfnhs.begin(), pfnhs.end()),
194  boost::make_filter_iterator(
195  isolationConeFilter,
196  pfnhs.end(), pfnhs.end())
197  );
198 
199  tau.addPiZeros(
201  boost::make_filter_iterator(
202  isolationConeFilterPiZero,
203  piZeroCombo->remainder_begin(), piZeroCombo->remainder_end()),
204  boost::make_filter_iterator(
205  isolationConeFilterPiZero,
206  piZeroCombo->remainder_end(), piZeroCombo->remainder_end())
207  );
208 
209  tau.addPiZeros(
211  boost::make_filter_iterator(
212  isolationConeFilterPiZero,
213  piZero_end, piZeros.end()),
214  boost::make_filter_iterator(
215  isolationConeFilterPiZero,
216  piZeros.end(), piZeros.end())
217  );
218 
219  output.push_back(tau.get(true));
220  }
221  }
222  }
223  return output.release();
224 }
225 }} // end namespace reco::tau
226 
230  "RecoTauBuilderCombinatoricPlugin");
T getParameter(std::string const &) const
const edm::Handle< PFCandidateCollection > & getPFCands() const
Hack to be able to convert Ptrs to Refs.
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
std::auto_ptr< reco::PFTau > get(bool setupLeadingCandidates=true)
void reserve(Region region, ParticleType type, size_t size)
Reserve a set amount of space for a given RefVector.
void addPiZeros(Region region, InputIterator begin, InputIterator end)
Add a list of pi zeros to the input collection.
std::vector< reco::PFCandidatePtr > PFCandPtrs
boost::ptr_vector< reco::PFTau > output_type
std::auto_ptr< output_type > return_type
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
virtual return_type operator()(const reco::PFJetRef &jet, const std::vector< RecoTauPiZero > &piZeros) const
Coll filterRefs(const Coll &refcoll) const
Filter a ref vector of PFCandidates.
tuple pset
Definition: CrabTask.py:85
const reco::Candidate::LorentzVector & p4() const
void reservePiZero(Region region, size_t size)
Reserve a set amount of space for the PiZeros.
std::vector< PFCandidatePtr > pfChargedCands(const PFJet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
void addPFCands(Region region, ParticleType type, InputIterator begin, InputIterator end)
#define DEFINE_EDM_PLUGIN(factory, type, name)
const reco::VertexRef & primaryVertex() const
Get primary vertex associated to this event.