CMS 3D CMS Logo

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  edm::InputTag srcWeights = iConfig.getParameter<edm::InputTag>("srcWeights");
22  if (!srcWeights.label().empty())
24 
25  for (std::vector<unsigned>::const_iterator n = Njets_.begin(); n != Njets_.end(); ++n) {
26  std::ostringstream ecfN_str;
27  std::shared_ptr<fastjet::FunctionOfPseudoJet<double>> pfunc;
28 
29  if (ecftype_ == "ECF" || ecftype_.empty()) {
30  ecfN_str << "ecf" << *n;
31  pfunc.reset(new fastjet::contrib::EnergyCorrelator(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
32  } else if (ecftype_ == "C") {
33  ecfN_str << "ecfC" << *n;
34  pfunc.reset(new fastjet::contrib::EnergyCorrelatorCseries(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
35  } else if (ecftype_ == "D") {
36  ecfN_str << "ecfD" << *n;
37  pfunc.reset(
38  new fastjet::contrib::EnergyCorrelatorGeneralizedD2(alpha_, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
39  } else if (ecftype_ == "N") {
40  ecfN_str << "ecfN" << *n;
41  pfunc.reset(new fastjet::contrib::EnergyCorrelatorNseries(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
42  } else if (ecftype_ == "M") {
43  ecfN_str << "ecfM" << *n;
44  pfunc.reset(new fastjet::contrib::EnergyCorrelatorMseries(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
45  } else if (ecftype_ == "U") {
46  ecfN_str << "ecfU" << *n;
47  pfunc.reset(new fastjet::contrib::EnergyCorrelatorUseries(*n, beta_, fastjet::contrib::EnergyCorrelator::pt_R));
48  }
49  variables_.push_back(ecfN_str.str());
50  produces<edm::ValueMap<float>>(ecfN_str.str());
51  routine_.push_back(pfunc);
52 
54  }
55 }
56 
58  // read input collection
60  iEvent.getByToken(src_token_, jets);
61 
62  // Get Weights Collection
65 
66  unsigned i = 0;
67  for (std::vector<unsigned>::const_iterator n = Njets_.begin(); n != Njets_.end(); ++n) {
68  // prepare room for output
69  std::vector<float> ecfN;
70  ecfN.reserve(jets->size());
71 
72  for (typename edm::View<reco::Jet>::const_iterator jetIt = jets->begin(); jetIt != jets->end(); ++jetIt) {
73  edm::Ptr<reco::Jet> jetPtr = jets->ptrAt(jetIt - jets->begin());
74 
75  float t = -1.0;
76  if (selectors_[n - Njets_.begin()](*jetIt))
77  t = getECF(i, jetPtr);
78 
79  ecfN.push_back(t);
80  }
81 
82  auto outT = std::make_unique<edm::ValueMap<float>>();
83  edm::ValueMap<float>::Filler fillerT(*outT);
84  fillerT.insert(jets, ecfN.begin(), ecfN.end());
85  fillerT.fill();
86 
87  iEvent.put(std::move(outT), variables_[i]);
88  ++i;
89  }
90 }
91 
92 float ECFAdder::getECF(unsigned index, const edm::Ptr<reco::Jet>& object) const {
93  std::vector<fastjet::PseudoJet> fjParticles;
94  for (unsigned k = 0; k < object->numberOfDaughters(); ++k) {
95  const reco::CandidatePtr& dp = object->daughterPtr(k);
96  if (dp.isNonnull() && dp.isAvailable()) {
97  // Here, the daughters are the "end" node, so this is a PFJet
98  if (dp->numberOfDaughters() == 0) {
99  if (object->isWeighted()) {
101  throw cms::Exception("MissingConstituentWeight")
102  << "ECFAdder: No weights (e.g. PUPPI) given for weighted jet collection" << std::endl;
103  float w = (*weightsHandle_)[dp];
104  fjParticles.push_back(fastjet::PseudoJet(dp->px() * w, dp->py() * w, dp->pz() * w, dp->energy() * w));
105  } else
106  fjParticles.push_back(fastjet::PseudoJet(dp->px(), dp->py(), dp->pz(), dp->energy()));
107  } else { // Otherwise, this is a BasicJet, so you need to descend further.
108  auto subjet = dynamic_cast<reco::Jet const*>(dp.get());
109  for (unsigned l = 0; l < subjet->numberOfDaughters(); ++l) {
110  if (subjet != nullptr) {
111  const reco::CandidatePtr& ddp = subjet->daughterPtr(l);
112  if (subjet->isWeighted()) {
114  throw cms::Exception("MissingConstituentWeight")
115  << "ECFAdder: No weights (e.g. PUPPI) given for weighted jet collection" << std::endl;
116  float w = (*weightsHandle_)[ddp];
117  fjParticles.push_back(fastjet::PseudoJet(ddp->px() * w, ddp->py() * w, ddp->pz() * w, ddp->energy() * w));
118  } else
119  fjParticles.push_back(fastjet::PseudoJet(ddp->px(), ddp->py(), ddp->pz(), ddp->energy()));
120  } else {
121  edm::LogWarning("MissingJetConstituent") << "BasicJet constituent required for ECF computation is missing!";
122  }
123  }
124  } // end if basic jet
125  } // end if daughter pointer is nonnull and available
126  else
127  edm::LogWarning("MissingJetConstituent") << "Jet constituent required for ECF computation is missing!";
128  }
129  if (fjParticles.size() > Njets_[index]) {
130  return routine_[index]->result(join(fjParticles));
131  } else {
132  return -1.0;
133  }
134 }
135 
136 // ParameterSet description for module
139  iDesc.setComment("Energy Correlation Functions adder");
140 
141  iDesc.add<edm::InputTag>("src", edm::InputTag("no default"))->setComment("input collection");
142  iDesc.add<std::vector<unsigned>>("Njets", {1, 2, 3})->setComment("Number of jets to emulate");
143  iDesc.add<std::vector<std::string>>("cuts", {"", "", ""})
144  ->setComment("Jet selections for each N value. Size must match Njets.");
145  iDesc.add<double>("alpha", 1.0)->setComment("alpha factor, only valid for N2");
146  iDesc.add<double>("beta", 1.0)->setComment("angularity factor");
147  iDesc.add<std::string>("ecftype", "")->setComment("ECF type: ECF or empty; C; D; N; M; U;");
148  iDesc.add<edm::InputTag>("srcWeights", edm::InputTag("puppi"));
149  descriptions.add("ECFAdder", iDesc);
150 }
151 
ConfigurationDescriptions.h
ECFAdder::beta_
double beta_
Definition: ECFAdder.h:32
resolutioncreator_cfi.object
object
Definition: resolutioncreator_cfi.py:4
mps_fire.i
i
Definition: mps_fire.py:428
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
sistrip::View
View
Definition: ConstantsForView.h:26
ECFAdder::weightsHandle_
edm::ValueMap< float > const * weightsHandle_
Definition: ECFAdder.h:37
ECFAdder::variables_
std::vector< std::string > variables_
Definition: ECFAdder.h:30
ECFAdder::src_token_
edm::EDGetTokenT< edm::View< reco::Jet > > src_token_
Definition: ECFAdder.h:26
ECFAdder::getECF
float getECF(unsigned index, const edm::Ptr< reco::Jet > &object) const
Definition: ECFAdder.cc:92
edm
HLT enums.
Definition: AlignableModifier.h:19
ECFAdder
Definition: ECFAdder.h:15
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
join
static std::string join(char **cmd)
Definition: RemoteFile.cc:17
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
edm::EDGetTokenT::isUninitialized
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:99
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:46
edm::Handle
Definition: AssociativeIterator.h:50
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
ECFAdder::cuts_
std::vector< std::string > cuts_
Definition: ECFAdder.h:28
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Jet
Definition: Jet.py:1
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
Calorimetry_cff.dp
dp
Definition: Calorimetry_cff.py:158
w
const double w
Definition: UKUtility.cc:23
ECFAdder.h
ECFAdder::routine_
std::vector< std::shared_ptr< fastjet::FunctionOfPseudoJet< double > > > routine_
Definition: ECFAdder.h:34
dqmdumpme.k
k
Definition: dqmdumpme.py:60
ParameterSetDescription.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
ECFAdder::ECFAdder
ECFAdder(const edm::ParameterSet &iConfig)
Definition: ECFAdder.cc:9
ECFAdder::selectors_
std::vector< StringCutObjectSelector< reco::Jet > > selectors_
Definition: ECFAdder.h:35
edm::ParameterSet
Definition: ParameterSet.h:47
edm::ParameterSetDescription::setComment
void setComment(std::string const &value)
Definition: ParameterSetDescription.cc:33
iEvent
int iEvent
Definition: GenABIO.cc:224
ECFAdder::produce
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Definition: ECFAdder.cc:57
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
edm::EventSetup
Definition: EventSetup.h:58
ECFAdder::ecftype_
std::string ecftype_
Definition: ECFAdder.h:29
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HLT_FULL_cff.srcWeights
srcWeights
Definition: HLT_FULL_cff.py:8709
cmsLHEtoEOSManager.l
l
Definition: cmsLHEtoEOSManager.py:204
edm::Ptr< reco::Jet >
ECFAdder::input_weights_token_
edm::EDGetTokenT< edm::ValueMap< float > > input_weights_token_
Definition: ECFAdder.h:36
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
edm::ValueMap< float >
StringCutObjectSelector< reco::Jet >
ECFAdder::alpha_
double alpha_
Definition: ECFAdder.h:31
Exception
Definition: hltDiff.cc:245
ECFAdder::Njets_
std::vector< unsigned > Njets_
Definition: ECFAdder.h:27
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Exception.h
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
edm::View::const_iterator
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
cms::Exception
Definition: Exception.h:70
edm::helper::Filler
Definition: ValueMap.h:22
edm::ParameterDescriptionNode::setComment
void setComment(std::string const &value)
Definition: ParameterDescriptionNode.cc:106
edm::Event
Definition: Event.h:73
submitPVValidationJobs.t
string t
Definition: submitPVValidationJobs.py:644
edm::InputTag
Definition: InputTag.h:15
ECFAdder::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: ECFAdder.cc:137