48 #include "fastjet/ClusterSequence.hh"
49 #include "fastjet/JetDefinition.hh"
50 #include "fastjet/PseudoJet.hh"
74 using JetDefPtr = std::shared_ptr<fastjet::JetDefinition>;
144 std::map<int, std::array<double, 4>>
eTop4_;
146 typedef std::pair<double, double>
EtaPhi;
151 : dropZeroTowers_(iConfig.getParameter<
bool>(
"dropZeroTowers")),
152 medianWindowWidth_(iConfig.getParameter<
int>(
"medianWindowWidth")),
153 minimumTowersFraction_(iConfig.getParameter<double>(
"minimumTowersFraction")),
154 nSigmaPU_(iConfig.getParameter<double>(
"nSigmaPU")),
155 puPtMin_(iConfig.getParameter<double>(
"puPtMin")),
156 radiusPU_(iConfig.getParameter<double>(
"radiusPU")),
157 rParam_(iConfig.getParameter<double>(
"rParam")),
158 towSigmaCut_(iConfig.getParameter<double>(
"towSigmaCut")),
162 produces<std::vector<double>>(
"mapEtaEdges");
163 produces<std::vector<double>>(
"mapToRho");
164 produces<std::vector<double>>(
"mapToRhoMedian");
165 produces<std::vector<double>>(
"mapToRhoExtra");
166 produces<std::vector<double>>(
"mapToRhoM");
167 produces<std::vector<int>>(
"mapToNTow");
168 produces<std::vector<double>>(
"mapToTowExcludePt");
169 produces<std::vector<double>>(
"mapToTowExcludePhi");
170 produces<std::vector<double>>(
"mapToTowExcludeEta");
182 inputs_.reserve(inputView.size());
183 for (
auto const&
input : inputView)
211 std::vector<fastjet::PseudoJet> orphanInput;
237 LogDebug(
"PileUpSubtractor") <<
"The subtractor setting up geometry...\n";
241 int ietaold = -10000;
246 for (
auto const& did : alldid) {
251 if (hid.
ieta() != ietaold) {
252 ietaold = hid.
ieta();
277 LogDebug(
"PileUpSubtractor") <<
"The subtractor calculating pedestals...\n";
279 std::map<int, double> emean2;
280 std::map<int, int> ntowers;
282 int ietaold = -10000;
285 for (
int vi = 0; vi <
nEtaTow_; ++vi) {
304 eTop4_[
i] = {{0., 0., 0., 0.}};
307 for (
auto const& input_object : coll) {
309 double original_Et = originalTower->
et();
310 int ieta0 = originalTower->
ieta();
312 if (original_Et >
eTop4_[ieta0][0]) {
316 eTop4_[ieta0][0] = original_Et;
317 }
else if (original_Et >
eTop4_[ieta0][1]) {
320 eTop4_[ieta0][1] = original_Et;
321 }
else if (original_Et >
eTop4_[ieta0][2]) {
323 eTop4_[ieta0][2] = original_Et;
324 }
else if (original_Et >
eTop4_[ieta0][3]) {
325 eTop4_[ieta0][3] = original_Et;
329 emean2[ieta0] = emean2[ieta0] + original_Et * original_Et;
330 if (ieta0 - ietaold != 0) {
347 double e2 = emean2[it];
354 LogDebug(
"PileUpSubtractor") <<
" ieta: " << it <<
" number of towers: " <<
nt <<
" e1: " <<
e1 <<
" e2: " << e2
361 double eee = e2 / (double)
nt -
e1 *
e1 / (
double)(
nt *
nt);
368 for (
unsigned int lI = 0; lI < numToCheck; ++lI) {
389 int sign = (it < 0) ? -1 : 1;
417 LogDebug(
"PileUpSubtractor") <<
" ieta: " << it <<
" Pedestals: " <<
emean_[it] <<
" " <<
esigma_[it] <<
"\n";
424 LogDebug(
"PileUpSubtractor") <<
"The subtractor subtracting pedestals...\n";
426 std::vector<fastjet::PseudoJet> newcoll;
427 for (
auto& input_object : coll) {
428 int index = input_object.user_index();
430 int it = itow->
ieta();
432 double original_Et = itow->
et();
434 float mScale = etnew / input_object.Et();
439 input_object.px() * mScale, input_object.py() * mScale, input_object.pz() * mScale, input_object.e() * mScale);
441 input_object.reset_momentum(towP4.px(), towP4.py(), towP4.pz(), towP4.energy());
442 input_object.set_user_index(
index);
445 newcoll.push_back(input_object);
453 LogDebug(
"PileUpSubtractor") <<
"The subtractor calculating orphan input...\n";
457 std::vector<int> jettowers;
458 std::map<std::pair<int, int>,
int> excludedTowers;
459 for (
auto const& pseudojetTMP :
fjJets_) {
460 EtaPhi jet_etaphi(pseudojetTMP.eta(), pseudojetTMP.phi());
466 double dr2 =
reco::deltaR2(im.eta, im.phi, jet_etaphi.first, jet_etaphi.second);
470 excludedTowers[std::pair(im.ieta, im.iphi)] = 1;
477 int ie = originalTower->
ieta();
478 int ip = originalTower->
iphi();
480 if (excludedTowers[std::pair<int, int>(ie, ip)]) {
481 jettowers.push_back(
index);
489 auto itjet =
std::find(jettowers.begin(), jettowers.end(),
index);
490 if (itjet == jettowers.end()) {
491 orphanInput.emplace_back(originalTower->
px(), originalTower->
py(), originalTower->
pz(), originalTower->
energy());
492 orphanInput.back().set_user_index(
index);
506 std::vector<std::pair<std::size_t, double>>
order;
507 for (std::size_t
i = 0;
i <
size; ++
i) {
512 order.begin(),
order.end(), [](
auto const& pair0,
auto const& pair1) {
return pair0.second < pair1.second; });
514 std::vector<double> sortedEtaEdgeLow(
size);
515 std::vector<double> sortedEtaEdgeHigh(
size);
517 auto mapToRhoOut = std::make_unique<std::vector<double>>(
size);
518 auto mapToRhoExtraOut = std::make_unique<std::vector<double>>(
size);
519 auto mapToRhoMOut = std::make_unique<std::vector<double>>(
size);
520 auto mapToNTowOut = std::make_unique<std::vector<int>>(
size);
522 for (std::size_t
i = 0;
i <
size; ++
i) {
534 auto mapToRhoMedianOut = std::make_unique<std::vector<double>>(
size);
537 auto centre = mapToRhoOut->begin() +
index;
539 std::nth_element(rhoRange.begin(), rhoRange.begin() +
medianWindowWidth_, rhoRange.end());
543 auto mapEtaRangesOut = std::make_unique<std::vector<double>>();
545 mapEtaRangesOut->assign(sortedEtaEdgeLow.begin(), sortedEtaEdgeLow.end());
546 mapEtaRangesOut->push_back(sortedEtaEdgeHigh.back());
567 desc.add<
int>(
"medianWindowWidth", 2);
568 desc.add<
double>(
"towSigmaCut", 5.0);
569 desc.add<
double>(
"puPtMin", 15.0);
570 desc.add<
double>(
"rParam", 0.3);
571 desc.add<
double>(
"nSigmaPU", 1.0);
572 desc.add<
double>(
"radiusPU", 0.5);
573 desc.add<
double>(
"minimumTowersFraction", 0.7);
574 desc.add<
bool>(
"dropZeroTowers",
true);
575 descriptions.
add(
"hiPuRhoProducer",
desc);