CMS 3D CMS Logo

BetaStarVarProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PhysicsTools/NanoAOD
4 // Class: BetaStarVarProducer
5 //
18 //
19 // Original Author: Sal Rappoccio
20 // Further editing: Hannu Siikonen
21 //
22 //
23 
24 #include <memory>
25 
28 
31 
34 
38 
40 
41 template <typename T>
43 public:
44  explicit BetaStarVarProducer(const edm::ParameterSet &iConfig)
45  : srcJet_(consumes<edm::View<pat::Jet>>(iConfig.getParameter<edm::InputTag>("srcJet"))),
46  srcPF_(consumes<edm::View<T>>(iConfig.getParameter<edm::InputTag>("srcPF"))),
47  maxDR_(iConfig.getParameter<double>("maxDR")) {
48  produces<edm::ValueMap<float>>("chargedFromPV0EnergyFraction");
49  produces<edm::ValueMap<float>>("chargedFromPV1EnergyFraction");
50  produces<edm::ValueMap<float>>("chargedFromPV2EnergyFraction");
51  produces<edm::ValueMap<float>>("chargedFromPV3EnergyFraction");
52  }
53  ~BetaStarVarProducer() override{};
54 
55  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
56 
57 private:
58  void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
59 
60  void fillValueMaps(std::string valueName,
61  const std::vector<float> &values,
63  edm::Event &iEvent) const;
64  std::tuple<float, float, float, float> calculateCHSEnergies(edm::Ptr<pat::Jet> const &jet, double chefrompv0) const;
65 
68  double maxDR_;
69 };
70 
71 template <typename T>
74  iEvent.getByToken(srcJet_, srcJet);
76  iEvent.getByToken(srcPF_, srcPF);
77 
78  unsigned int nJet = srcJet->size();
79 
80  // Create an injective mapping from jets to charged PF candidates with fromPV==0 within the jet cone.
81  // These are the PF candidates supposedly removed by the CHS method - i.e. removed charged pileup.
82  std::vector<double> jet2pue(nJet, 0);
83  for (const auto &cand : *srcPF) {
84  if (cand.charge() == 0 or cand.fromPV() != 0)
85  continue;
86  // Looking for the closest match for this charged frompv==0 candidate
87  float dRMin = 999.;
88  int bestjet = -1, jetidx = 0;
89  for (const auto &jet : *srcJet) {
90  float dR = reco::deltaR(jet, cand);
91  if (dR < dRMin) {
92  dRMin = dR;
93  bestjet = jetidx;
94  }
95  ++jetidx;
96  }
97  // If the candidate is closer than the jet radius to the jet axis, this is a PU particle of interest for the selected jet
98  if (dRMin < maxDR_)
99  jet2pue[bestjet] += cand.energy();
100  }
101 
102  std::vector<float> chargedFromPV0EnergyFraction(nJet, -1);
103  std::vector<float> chargedFromPV1EnergyFraction(nJet, -1);
104  std::vector<float> chargedFromPV2EnergyFraction(nJet, -1);
105  std::vector<float> chargedFromPV3EnergyFraction(nJet, -1);
106 
107  for (unsigned int ij = 0; ij < nJet; ++ij) {
108  auto vals = calculateCHSEnergies(srcJet->ptrAt(ij), jet2pue[ij]);
109  chargedFromPV0EnergyFraction[ij] = std::get<0>(vals);
110  chargedFromPV1EnergyFraction[ij] = std::get<1>(vals);
111  chargedFromPV2EnergyFraction[ij] = std::get<2>(vals);
112  chargedFromPV3EnergyFraction[ij] = std::get<3>(vals);
113  }
114 
115  fillValueMaps("chargedFromPV0EnergyFraction", chargedFromPV0EnergyFraction, srcJet, iEvent);
116  fillValueMaps("chargedFromPV1EnergyFraction", chargedFromPV1EnergyFraction, srcJet, iEvent);
117  fillValueMaps("chargedFromPV2EnergyFraction", chargedFromPV2EnergyFraction, srcJet, iEvent);
118  fillValueMaps("chargedFromPV3EnergyFraction", chargedFromPV3EnergyFraction, srcJet, iEvent);
119 }
120 
121 template <typename T>
123  const std::vector<float> &values,
125  edm::Event &iEvent) const {
126  std::unique_ptr<edm::ValueMap<float>> valuesPtr(new edm::ValueMap<float>());
127  edm::ValueMap<float>::Filler valuesFiller(*valuesPtr);
128  valuesFiller.insert(srcJet, values.begin(), values.end());
129  valuesFiller.fill();
130  iEvent.put(std::move(valuesPtr), valueName);
131 }
132 
133 template <typename T>
134 std::tuple<float, float, float, float> BetaStarVarProducer<T>::calculateCHSEnergies(edm::Ptr<pat::Jet> const &ijet,
135  double chefrompv0) const {
136  // Keep track of energy (present in the jet) for PU stuff
137  double chefrompv1 = 0.0;
138  double chefrompv2 = 0.0;
139  double chefrompv3 = 0.0;
140 
141  // Loop through the PF candidates within the jet.
142  // Store the sum of their energy, and their indices.
143  for (auto pidx = 0u; pidx < ijet->numberOfDaughters(); ++pidx) {
144  auto *dtr = dynamic_cast<const pat::PackedCandidate *>(ijet->daughter(pidx));
145  if (dtr->charge() == 0)
146  continue;
147  if (dtr->fromPV() == 1)
148  chefrompv1 += dtr->energy();
149  else if (dtr->fromPV() == 2)
150  chefrompv2 += dtr->energy();
151  else if (dtr->fromPV() == 3)
152  chefrompv3 += dtr->energy();
153  }
154 
155  // Now get the fractions relative to the raw jet.
156  auto rawP4 = ijet->correctedP4(0);
157  double chffpv0 = chefrompv0 / rawP4.energy();
158  double chffpv1 = chefrompv1 / rawP4.energy();
159  double chffpv2 = chefrompv2 / rawP4.energy();
160  double chffpv3 = chefrompv3 / rawP4.energy();
161 
162  return std::tuple<float, float, float, float>(chffpv0, chffpv1, chffpv2, chffpv3);
163 }
164 
165 template <typename T>
168  desc.add<edm::InputTag>("srcJet")->setComment("jet input collection");
169  desc.add<edm::InputTag>("srcPF")->setComment("PF candidate input collection");
170  desc.add<double>("maxDR")->setComment("Maximum DR to consider for jet-->pf cand association");
171  std::string modname("BetaStarVarProducer");
172  descriptions.add(modname, desc);
173 }
174 
176 
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
BetaStarVarProducer< pat::PackedCandidate > BetaStarPackedCandidateVarProducer
void fillValueMaps(std::string valueName, const std::vector< float > &values, edm::Handle< edm::View< pat::Jet >> &srcJet, edm::Event &iEvent) const
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition: HeavyIon.h:7
edm::EDGetTokenT< edm::View< pat::Jet > > srcJet_
size_t numberOfDaughters() const override
int iEvent
Definition: GenABIO.cc:224
Definition: Jet.py:1
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
edm::EDGetTokenT< edm::View< T > > srcPF_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
std::tuple< float, float, float, float > calculateCHSEnergies(edm::Ptr< pat::Jet > const &jet, double chefrompv0) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
BetaStarVarProducer(const edm::ParameterSet &iConfig)
HLT enums.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double energy() const override
energy
const reco::Candidate * daughter(size_t i) const override
long double T
const LorentzVector correctedP4(const std::string &level, const std::string &flavor="none", const std::string &set="") const
Definition: Jet.h:164
def move(src, dest)
Definition: eostools.py:511