CMS 3D CMS Logo

Phase1L1TJetCalibrator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1CaloTrigger
4 // Class: Phase1L1TJetCalibrator
5 //
17 //
18 // Original Author: Simone Bologna
19 // Created: Wed, 19 Dec 2018 12:44:23 GMT
20 //
21 // Rewrite: Vladimir Rekovic, Oct 2020
22 //
23 
24 // system include files
25 #include <memory>
26 
27 #include <algorithm>
28 #include <cmath>
29 
30 // user include files
33 
39 
42 
44 public:
46  ~Phase1L1TJetCalibrator() override;
47 
48  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
49 
50 private:
51  void produce(edm::Event&, const edm::EventSetup&) override;
52 
54  std::vector<double> absEtaBinning_;
55  size_t nBinsEta_;
56  std::vector<edm::ParameterSet> calibration_;
58 
59  std::vector<std::vector<double>> jetCalibrationFactorsBinnedInEta_;
60  std::vector<std::vector<double>> _jetCalibrationFactorsPtBins;
61 
62  // ----------member data ---------------------------
63 };
64 
65 //
66 // constants, enums and typedefs
67 //
68 
69 //
70 // static data member definitions
71 //
72 
73 //
74 // constructors and destructor
75 //
77  : inputCollectionTag_(edm::EDGetTokenT<std::vector<reco::CaloJet>>(
78  consumes<std::vector<reco::CaloJet>>(iConfig.getParameter<edm::InputTag>("inputCollectionTag")))),
79  absEtaBinning_(iConfig.getParameter<std::vector<double>>("absEtaBinning")),
80  nBinsEta_(absEtaBinning_.size() - 1),
81  calibration_(iConfig.getParameter<std::vector<edm::ParameterSet>>("calibration")),
82  outputCollectionName_(iConfig.getParameter<std::string>("outputCollectionName"))
83 
84 {
85  for (const auto& pset : calibration_) {
86  _jetCalibrationFactorsPtBins.emplace_back(pset.getParameter<std::vector<double>>("l1tPtBins"));
87  jetCalibrationFactorsBinnedInEta_.emplace_back(pset.getParameter<std::vector<double>>("l1tCalibrationFactors"));
88  }
89 
90  produces<std::vector<reco::CaloJet>>(outputCollectionName_).setBranchAlias(outputCollectionName_);
91 }
92 
94 
95 //
96 // member functions
97 //
98 
99 // ------------ method called to produce the data ------------
101  edm::Handle<std::vector<reco::CaloJet>> inputCollectionHandle;
102  iEvent.getByToken(inputCollectionTag_, inputCollectionHandle);
103 
104  auto calibratedCollectionPtr = std::make_unique<std::vector<reco::CaloJet>>();
105 
106  // for each candidate:
107  // 1 get pt and eta
108  // 2 run a dicotomic search on the eta vector to find the eta index
109  // 3 fetch the corresponding calibration elements
110  // 4 run a dicotomic search on the pt binning to find the pt index
111  // 5 fetch the calibration factor
112  // 6 update the candidate pt by applying the calibration factor
113  // 7 store calibrated candidate in a new collection
114 
115  calibratedCollectionPtr->reserve(inputCollectionHandle->size());
116 
117  for (const auto& candidate : *inputCollectionHandle) {
118  // 1
119  float pt = candidate.pt();
120  float eta = candidate.eta();
121 
122  //2
123  auto etaBin = upper_bound(absEtaBinning_.begin(), absEtaBinning_.end(), fabs(eta));
124  int etaIndex = etaBin - absEtaBinning_.begin() - 1;
125 
126  //3
127  const std::vector<double>& l1tPtBins = _jetCalibrationFactorsPtBins[etaIndex];
128  const std::vector<double>& l1tCalibrationFactors = jetCalibrationFactorsBinnedInEta_[etaIndex];
129 
130  //4
131  auto ptBin = upper_bound(l1tPtBins.begin(), l1tPtBins.end(), pt);
132  int ptBinIndex = ptBin - l1tPtBins.begin() - 1;
133 
134  //5
135  const double& l1tCalibrationFactor = l1tCalibrationFactors[ptBinIndex];
136 
137  //6 and 7
138  reco::Candidate::PolarLorentzVector candidateP4(candidate.polarP4());
139  candidateP4.SetPt(candidateP4.pt() * l1tCalibrationFactor);
140  calibratedCollectionPtr->emplace_back(candidate);
141  calibratedCollectionPtr->back().setP4(candidateP4);
142 
143 #ifdef DEBUG
144  if (newCandidate->pt() < 0) {
145  LogDebug("Phase1L1TJetCalibrator") << "######################" << std::endl;
146  LogDebug("Phase1L1TJetCalibrator") << "PRE-CALIBRATION " << std::endl;
147  LogDebug("Phase1L1TJetCalibrator") << "\t Jet properties (pt, eta, phi, pile-up): " << candidate.pt() << "\t"
148  << candidate.eta() << "\t" LogDebug("Phase1L1TJetCalibrator")
149  << candidate.phi() << "\t" << candidate.pileup() << std::endl;
150  LogDebug("Phase1L1TJetCalibrator") << "CALIBRATION " << std::endl;
151  LogDebug("Phase1L1TJetCalibrator") << "\t Using eta - pt - factor " << *etaBin << " - " << *ptBin << " - "
152  << l1tCalibrationFactor LogDebug("Phase1L1TJetCalibrator") << std::endl;
153  LogDebug("Phase1L1TJetCalibrator") << "POST-CALIBRATION " << std::endl;
154  LogDebug("Phase1L1TJetCalibrator") << "\t Jet properties (pt, eta, phi, pile-up): " << newCandidate->pt() << "\t"
155  << newCandidate->eta() << "\t" << newCandidate->phi() << "\t"
156  << newCandidate->pileup() << std::endl;
157  }
158 #endif
159  }
160 
161  // finally, sort the collection by pt
162  std::sort(calibratedCollectionPtr->begin(),
163  calibratedCollectionPtr->end(),
164  [](const reco::CaloJet& jet1, const reco::CaloJet& jet2) { return jet1.pt() > jet2.pt(); });
165 
166  iEvent.put(std::move(calibratedCollectionPtr), outputCollectionName_);
167 }
168 
171  desc.add<edm::InputTag>("inputCollectionTag",
172  edm::InputTag("l1tPhase1JetProducer", "UncalibratedPhase1L1TJetFromPfCandidates"));
173  desc.add<std::vector<double>>("absEtaBinning");
174  std::vector<edm::ParameterSet> vDefaults;
176  validator.add<double>("etaMax");
177  validator.add<double>("etaMin");
178  validator.add<std::vector<double>>("l1tCalibrationFactors");
179  validator.add<std::vector<double>>("l1tPtBins");
180  desc.addVPSet("calibration", validator, vDefaults);
181  desc.add<std::string>("outputCollectionName", "Phase1L1TJetFromPfCandidates");
182  descriptions.add("Phase1L1TJetCalibrator", desc);
183 }
184 
185 //define this as a plug-in
size
Write out results.
Jets made from CaloTowers.
Definition: CaloJet.h:27
void produce(edm::Event &, const edm::EventSetup &) override
double pt() const final
transverse momentum
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< edm::ParameterSet > calibration_
int iEvent
Definition: GenABIO.cc:224
std::vector< std::vector< double > > _jetCalibrationFactorsPtBins
Phase1L1TJetCalibrator(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< double > absEtaBinning_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
std::vector< std::vector< double > > jetCalibrationFactorsBinnedInEta_
edm::EDGetTokenT< std::vector< reco::CaloJet > > inputCollectionTag_
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38