CMS 3D CMS Logo

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  */
11 
12 #include <algorithm>
13 
19 
23 
25 
26 namespace reco { namespace tau {
27 
29  public:
32  // Return type is auto_ptr<PiZeroVector>
33  return_type operator()(const reco::PFJet& jet) const override;
34 
35  private:
37  double minMass_;
38  double maxMass_;
39  unsigned int maxInputGammas_;
40  unsigned int choose_;
42 };
43 
47  "qualityCuts").getParameterSet("signalQualityCuts")) {
48  minMass_ = pset.getParameter<double>("minMass");
49  maxMass_ = pset.getParameter<double>("maxMass");
50  maxInputGammas_ = pset.getParameter<unsigned int>("maxInputGammas");
51  choose_ = pset.getParameter<unsigned int>("choose");
52 }
53 
56  const reco::PFJet& jet) const {
57  // Get list of gamma candidates
58  typedef std::vector<reco::PFCandidatePtr> PFCandPtrs;
59  typedef PFCandPtrs::const_iterator PFCandIter;
61 
62  PFCandPtrs pfGammaCands = qcuts_.filterCandRefs(pfGammas(jet));
63  // Check if we have anything to do...
64  if (pfGammaCands.size() < choose_)
65  return output.release();
66 
67  // Define the valid range of gammas to use
68  PFCandIter start_iter = pfGammaCands.begin();
69  PFCandIter end_iter = pfGammaCands.end();
70 
71  // Only take the desired number of piZeros
72  end_iter = takeNElements(start_iter, end_iter, maxInputGammas_);
73 
74  // Build the combinatoric generator
75  typedef CombinatoricGenerator<PFCandPtrs> ComboGenerator;
76  ComboGenerator generator(start_iter, end_iter, choose_);
77 
78  // Find all possible combinations
79  for (ComboGenerator::iterator combo = generator.begin();
80  combo != generator.end(); ++combo) {
81  const Candidate::LorentzVector totalP4;
82  std::auto_ptr<RecoTauPiZero> piZero(
83  new RecoTauPiZero(0, totalP4, Candidate::Point(0, 0, 0),
84  111, 10001, true, RecoTauPiZero::kCombinatoric));
85  // Add our daughters from this combination
86  for (ComboGenerator::combo_iterator candidate = combo->combo_begin();
87  candidate != combo->combo_end(); ++candidate) {
88  piZero->addDaughter(*candidate);
89  }
90  p4Builder_.set(*piZero);
91 
92  if (piZero->daughterPtr(0).isNonnull())
93  piZero->setVertex(piZero->daughterPtr(0)->vertex());
94 
95  if ((maxMass_ < 0 || piZero->mass() < maxMass_) &&
96  piZero->mass() > minMass_)
97  output.push_back(piZero);
98  }
99  return output.release();
100 }
101 
102 }} // end namespace reco::tau
103 
107  "RecoTauPiZeroCombinatoricPlugin");
T getParameter(std::string const &) const
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
RecoTauPiZeroCombinatoricPlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
Coll filterCandRefs(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:21
return_type operator()(const reco::PFJet &jet) const override
Build a collection of piZeros from objects in the input jet.
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
std::auto_ptr< PiZeroVector > return_type
fixed size matrix
boost::ptr_vector< RecoTauPiZero > PiZeroVector
PFCandPtrs::iterator PFCandIter
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
#define DEFINE_EDM_PLUGIN(factory, type, name)
void set(reco::Candidate &c) const
set up a candidate
def move(src, dest)
Definition: eostools.py:510