CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ECFAdder.cc
Go to the documentation of this file.
2 #include "fastjet/PseudoJet.hh"
3 #include "fastjet/ClusterSequence.hh"
8 
10  : src_(iConfig.getParameter<edm::InputTag>("src")),
11  src_token_(consumes<edm::View<reco::Jet>>(src_)),
12  Njets_(iConfig.getParameter<std::vector<unsigned>>("Njets")),
13  cuts_(iConfig.getParameter<std::vector<std::string>>("cuts")),
14  ecftype_(iConfig.getParameter<std::string>("ecftype")),
15  alpha_(iConfig.getParameter<double>("alpha")),
16  beta_(iConfig.getParameter<double>("beta")) {
17  if (cuts_.size() != Njets_.size()) {
18  throw cms::Exception("ConfigurationError") << "cuts and Njets must be the same size in ECFAdder" << std::endl;
19  }
20 
21  for (std::vector<unsigned>::const_iterator n = Njets_.begin(); n != Njets_.end(); ++n) {
22  std::ostringstream ecfN_str;
23  std::shared_ptr<fastjet::FunctionOfPseudoJet<double>> pfunc;
24 
25  if (ecftype_ == "ECF" || ecftype_.empty()) {
26  ecfN_str << "ecf" << *n;
27  pfunc.reset(new fastjet::contrib::EnergyCorrelator(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
28  } else if (ecftype_ == "C") {
29  ecfN_str << "ecfC" << *n;
30  pfunc.reset(new fastjet::contrib::EnergyCorrelatorCseries(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
31  } else if (ecftype_ == "D") {
32  ecfN_str << "ecfD" << *n;
33  pfunc.reset(
34  new fastjet::contrib::EnergyCorrelatorGeneralizedD2(alpha_, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
35  } else if (ecftype_ == "N") {
36  ecfN_str << "ecfN" << *n;
37  pfunc.reset(new fastjet::contrib::EnergyCorrelatorNseries(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
38  } else if (ecftype_ == "M") {
39  ecfN_str << "ecfM" << *n;
40  pfunc.reset(new fastjet::contrib::EnergyCorrelatorMseries(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
41  } else if (ecftype_ == "U") {
42  ecfN_str << "ecfU" << *n;
43  pfunc.reset(new fastjet::contrib::EnergyCorrelatorUseries(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
44  }
45  variables_.push_back(ecfN_str.str());
46  produces<edm::ValueMap<float>>(ecfN_str.str());
47  routine_.push_back(pfunc);
48 
50  }
51 }
52 
54  // read input collection
56  iEvent.getByToken(src_token_, jets);
57 
58  unsigned i = 0;
59  for (std::vector<unsigned>::const_iterator n = Njets_.begin(); n != Njets_.end(); ++n) {
60  // prepare room for output
61  std::vector<float> ecfN;
62  ecfN.reserve(jets->size());
63 
64  for (typename edm::View<reco::Jet>::const_iterator jetIt = jets->begin(); jetIt != jets->end(); ++jetIt) {
65  edm::Ptr<reco::Jet> jetPtr = jets->ptrAt(jetIt - jets->begin());
66 
67  float t = -1.0;
68  if (selectors_[n - Njets_.begin()](*jetIt))
69  t = getECF(i, jetPtr);
70 
71  ecfN.push_back(t);
72  }
73 
74  auto outT = std::make_unique<edm::ValueMap<float>>();
75  edm::ValueMap<float>::Filler fillerT(*outT);
76  fillerT.insert(jets, ecfN.begin(), ecfN.end());
77  fillerT.fill();
78 
79  iEvent.put(std::move(outT), variables_[i]);
80  ++i;
81  }
82 }
83 
84 float ECFAdder::getECF(unsigned index, const edm::Ptr<reco::Jet>& object) const {
85  std::vector<fastjet::PseudoJet> FJparticles;
86  for (unsigned k = 0; k < object->numberOfDaughters(); ++k) {
87  const reco::CandidatePtr& dp = object->daughterPtr(k);
88  if (dp.isNonnull() && dp.isAvailable()) {
89  // Here, the daughters are the "end" node, so this is a PFJet
90  if (dp->numberOfDaughters() == 0) {
91  FJparticles.push_back(fastjet::PseudoJet(dp->px(), dp->py(), dp->pz(), dp->energy()));
92  } else { // Otherwise, this is a BasicJet, so you need to descend further.
93  auto subjet = dynamic_cast<reco::Jet const*>(dp.get());
94  for (unsigned l = 0; l < subjet->numberOfDaughters(); ++l) {
95  if (subjet != nullptr) {
96  const reco::CandidatePtr& ddp = subjet->daughterPtr(l);
97  FJparticles.push_back(fastjet::PseudoJet(ddp->px(), ddp->py(), ddp->pz(), ddp->energy()));
98  } else {
99  edm::LogWarning("MissingJetConstituent") << "BasicJet constituent required for ECF computation is missing!";
100  }
101  }
102  } // end if basic jet
103  } // end if daughter pointer is nonnull and available
104  else
105  edm::LogWarning("MissingJetConstituent") << "Jet constituent required for ECF computation is missing!";
106  }
107  if (FJparticles.size() > Njets_[index]) {
108  return routine_[index]->result(join(FJparticles));
109  } else {
110  return -1.0;
111  }
112 }
113 
114 // ParameterSet description for module
117  iDesc.setComment("Energy Correlation Functions adder");
118 
119  iDesc.add<edm::InputTag>("src", edm::InputTag("no default"))->setComment("input collection");
120  iDesc.add<std::vector<unsigned>>("Njets", {1, 2, 3})->setComment("Number of jets to emulate");
121  iDesc.add<std::vector<std::string>>("cuts", {"", "", ""})
122  ->setComment("Jet selections for each N value. Size must match Njets.");
123  iDesc.add<double>("alpha", 1.0)->setComment("alpha factor, only valid for N2");
124  iDesc.add<double>("beta", 1.0)->setComment("angularity factor");
125  iDesc.add<std::string>("ecftype", "")->setComment("ECF type: ECF or empty; C; D; N; M; U;");
126  descriptions.add("ECFAdder", iDesc);
127 }
128 
void setComment(std::string const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: ECFAdder.cc:115
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
ECFAdder(const edm::ParameterSet &iConfig)
Definition: ECFAdder.cc:9
std::string ecftype_
Definition: ECFAdder.h:29
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:139
Base class for all types of Jets.
Definition: Jet.h:20
bool isAvailable() const
Definition: Ptr.h:229
float getECF(unsigned index, const edm::Ptr< reco::Jet > &object) const
Definition: ECFAdder.cc:84
std::vector< std::string > variables_
Definition: ECFAdder.h:30
std::vector< std::string > cuts_
Definition: ECFAdder.h:28
std::vector< std::shared_ptr< fastjet::FunctionOfPseudoJet< double > > > routine_
Definition: ECFAdder.h:34
void setComment(std::string const &value)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double beta_
Definition: ECFAdder.h:32
Definition: Jet.py:1
edm::EDGetTokenT< edm::View< reco::Jet > > src_token_
Definition: ECFAdder.h:26
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< unsigned > Njets_
Definition: ECFAdder.h:27
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
static std::string join(char **cmd)
Definition: RemoteFile.cc:17
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Definition: ECFAdder.cc:53
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
double alpha_
Definition: ECFAdder.h:31
std::vector< StringCutObjectSelector< reco::Jet > > selectors_
Definition: ECFAdder.h:35
def move(src, dest)
Definition: eostools.py:511