CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauPiZeroCombinatoricPlugin.cc
Go to the documentation of this file.
1 /*
2  * RecoTauPiZeroCombinatoricPlugin
3  *
4  * Author: Evan K. Friis, UC Davis
5  *
6  * Build PiZero candidates out of all possible sets of <choose> gammas that are
7  * contained in the input PFJet. Optionally, the pi zero candidates are
8  * filtered by a min and max selection on their invariant mass.
9  *
10  * $Id $
11  */
12 
13 #include <algorithm>
14 
20 
24 
26 
27 namespace reco { namespace tau {
28 
30  public:
33  // Return type is auto_ptr<PiZeroVector>
34  return_type operator()(const reco::PFJet& jet) const;
35 
36  private:
38  double minMass_;
39  double maxMass_;
40  unsigned int maxInputGammas_;
41  unsigned int choose_;
43 };
44 
47  qcuts_(pset.getParameterSet(
48  "qualityCuts").getParameterSet("signalQualityCuts")) {
49  minMass_ = pset.getParameter<double>("minMass");
50  maxMass_ = pset.getParameter<double>("maxMass");
51  maxInputGammas_ = pset.getParameter<unsigned int>("maxInputGammas");
52  choose_ = pset.getParameter<unsigned int>("choose");
53 }
54 
57  const reco::PFJet& jet) const {
58  // Get list of gamma candidates
59  typedef std::vector<reco::PFCandidatePtr> PFCandPtrs;
60  typedef PFCandPtrs::const_iterator PFCandIter;
62 
63  PFCandPtrs pfGammaCands = qcuts_.filterRefs(pfGammas(jet));
64  // Check if we have anything to do...
65  if (pfGammaCands.size() < choose_)
66  return output.release();
67 
68  // Define the valid range of gammas to use
69  PFCandIter start_iter = pfGammaCands.begin();
70  PFCandIter end_iter = pfGammaCands.end();
71 
72  // Only take the desired number of piZeros
73  end_iter = takeNElements(start_iter, end_iter, maxInputGammas_);
74 
75  // Build the combinatoric generator
76  typedef CombinatoricGenerator<PFCandPtrs> ComboGenerator;
77  ComboGenerator generator(start_iter, end_iter, choose_);
78 
79  // Find all possible combinations
80  for (ComboGenerator::iterator combo = generator.begin();
81  combo != generator.end(); ++combo) {
82  const Candidate::LorentzVector totalP4;
83  std::auto_ptr<RecoTauPiZero> piZero(
84  new RecoTauPiZero(0, totalP4, Candidate::Point(0, 0, 0),
85  111, 10001, true, RecoTauPiZero::kCombinatoric));
86  // Add our daughters from this combination
87  for (ComboGenerator::combo_iterator candidate = combo->combo_begin();
88  candidate != combo->combo_end(); ++candidate) {
89  piZero->addDaughter(*candidate);
90  }
91  p4Builder_.set(*piZero);
92 
93  if (piZero->daughterPtr(0).isNonnull())
94  piZero->setVertex(piZero->daughterPtr(0)->vertex());
95 
96  if ((maxMass_ < 0 || piZero->mass() < maxMass_) &&
97  piZero->mass() > minMass_)
98  output.push_back(piZero);
99  }
100  return output.release();
101 }
102 
103 }} // end namespace reco::tau
104 
108  "RecoTauPiZeroCombinatoricPlugin");
T getParameter(std::string const &) const
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
Coll filterRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of PFCandidates.
std::vector< PFCandidatePtr > pfGammas(const PFJet &jet, bool sort=true)
Extract all pfGammas from a PFJet.
ParameterSet const & getParameterSet(ParameterSetID const &id)
std::vector< reco::PFCandidatePtr > PFCandPtrs
Jets made from PFObjects.
Definition: PFJet.h:22
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:38
std::auto_ptr< PiZeroVector > return_type
boost::ptr_vector< RecoTauPiZero > PiZeroVector
tuple mass
Definition: scaleCards.py:27
return_type operator()(const reco::PFJet &jet) const
Build a collection of piZeros from objects in the input jet.
PFCandPtrs::iterator PFCandIter
math::XYZPoint Point
point in the space
Definition: Candidate.h:42
#define DEFINE_EDM_PLUGIN(factory, type, name)
void set(reco::Candidate &c) const
set up a candidate