|
![CMS Logo](/cmsdoxygen/common/rightImage.jpg) |
#include <RecoHI/HiJetAlgos/plugins/HiPuRhoProducer.cc>
|
typedef std::pair< double, double > | EtaPhi |
|
|
std::vector< HcalDetId > | allgeomid_ |
|
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > | caloGeometryToken_ |
|
const edm::EDGetTokenT< CaloTowerCollection > | caloTowerToken_ |
|
const bool | dropZeroTowers_ |
|
std::map< int, double > | emean_ |
|
std::map< int, double > | esigma_ |
|
std::vector< double > | etaEdgeHi_ |
|
std::vector< double > | etaEdgeLow_ |
|
std::vector< double > | etaEdges_ |
|
std::map< int, std::array< double, 4 > > | eTop4_ |
|
ClusterSequencePtr | fjClusterSeq_ |
|
std::vector< fastjet::PseudoJet > | fjInputs_ |
|
JetDefPtr | fjJetDefinition_ |
|
std::vector< fastjet::PseudoJet > | fjJets_ |
|
std::vector< fastjet::PseudoJet > | fjOriginalInputs_ |
|
CaloGeometry const * | geo_ = nullptr |
|
std::map< int, int > | geomtowers_ |
|
int | ietamax_ |
|
int | ietamin_ |
|
const int | initialValue_ = -99 |
|
std::vector< const CaloTower * > | inputs_ |
|
const int | medianWindowWidth_ |
|
const double | minimumTowersFraction_ |
|
const double | nSigmaPU_ |
|
std::vector< int > | nTow_ |
|
std::map< int, int > | ntowersWithJets_ |
|
bool | postOrphan_ |
|
const double | puPtMin_ |
|
const double | radiusPU_ |
|
std::vector< double > | rho_ |
|
std::vector< double > | rhoExtra_ |
|
std::vector< double > | rhoM_ |
|
const double | rParam_ |
|
bool | setInitialValue_ |
|
std::vector< EtaPhiTower > | towermap_ |
|
std::vector< double > | towExcludeEta_ |
|
std::vector< double > | towExcludePhi_ |
|
std::vector< double > | towExcludePt_ |
|
const double | towSigmaCut_ |
|
std::array< float, nEtaTow_ > | vmean0_ |
|
std::array< float, nEtaTow_ > | vmean1_ |
|
std::array< int, nEtaTow_ > | vngeom_ |
|
std::array< int, nEtaTow_ > | vntow_ |
|
std::array< float, nEtaTow_ > | vrho0_ |
|
std::array< float, nEtaTow_ > | vrho1_ |
|
std::array< float, nEtaTow_ > | vrms0_ |
|
std::array< float, nEtaTow_ > | vrms1_ |
|
Description: Producer to dump Pu-jet style rho into event content Implementation: Just see MultipleAlgoIterator - re-implenting for use in CS jet with sigma subtraction and zeroing
Definition at line 69 of file HiPuRhoProducer.cc.
◆ ClusterSequencePtr
◆ EtaPhi
◆ JetDefPtr
◆ HiPuRhoProducer()
Definition at line 150 of file HiPuRhoProducer.cc.
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");
◆ calculateOrphanInput()
void HiPuRhoProducer::calculateOrphanInput |
( |
std::vector< fastjet::PseudoJet > & |
orphanInput | ) |
|
|
virtual |
Definition at line 452 of file HiPuRhoProducer.cc.
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);
References reco::deltaR2(), reco::LeafCandidate::energy(), reco::LeafCandidate::eta(), spr::find(), fjInputs_, fjJets_, fjOriginalInputs_, CaloTower::ieta(), input, inputs_, CaloTower::iphi(), LogDebug, minimumTowersFraction_, ntowersWithJets_, reco::LeafCandidate::phi(), postOrphan_, reco::LeafCandidate::pt(), puPtMin_, reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), radiusPU_, towermap_, towExcludeEta_, towExcludePhi_, and towExcludePt_.
Referenced by produce().
◆ calculatePedestal()
void HiPuRhoProducer::calculatePedestal |
( |
std::vector< fastjet::PseudoJet > const & |
coll | ) |
|
|
virtual |
Definition at line 276 of file HiPuRhoProducer.cc.
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";
References funct::abs(), StorageManager_cfg::e1, emean_, esigma_, CaloTower::et(), hi::etaedge, etaEdgeHi_, etaEdgeLow_, eTop4_, geomtowers_, submitPVValidationJobs::gt, mps_fire::i, CaloTower::ieta(), ietamax_, ietamin_, initialValue_, inputs_, LogDebug, M_PI, min(), minimumTowersFraction_, nEtaTow_, nSigmaPU_, nt, nTow_, ntowersWithJets_, postOrphan_, rho_, rhoExtra_, rhoM_, setInitialValue_, Validation_hcalonly_cfi::sign, findQualityFiles::size, mathSSE::sqrt(), towSigmaCut_, vmean0_, vmean1_, vngeom_, vntow_, vrho0_, vrho1_, vrms0_, and vrms1_.
Referenced by produce().
◆ fillDescriptions()
Definition at line 564 of file HiPuRhoProducer.cc.
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);
References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, and HLT_FULL_cff::InputTag.
◆ inputTowers()
void HiPuRhoProducer::inputTowers |
( |
| ) |
|
|
privatevirtual |
◆ produce()
Definition at line 174 of file HiPuRhoProducer.cc.
182 inputs_.reserve(inputView.size());
183 for (
auto const&
input : inputView)
211 std::vector<fastjet::PseudoJet> orphanInput;
References calculateOrphanInput(), calculatePedestal(), caloTowerToken_, etaEdgeHi_, etaEdgeLow_, etaEdges_, fjClusterSeq_, fjInputs_, fjJetDefinition_, fjJets_, fjOriginalInputs_, mps_fire::i, ietamax_, ietamin_, iEvent, input, inputs_, inputTowers(), nTow_, ntowersWithJets_, puPtMin_, putRho(), rho_, rhoExtra_, rhoM_, rParam_, setInitialValue_, setupGeometryMap(), subtractPedestal(), towExcludeEta_, towExcludePhi_, and towExcludePt_.
◆ putRho()
Definition at line 503 of file HiPuRhoProducer.cc.
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());
References etaEdgeHi_, etaEdgeLow_, mps_fire::i, iEvent, medianWindowWidth_, eostools::move(), nTow_, eventshapeDQM_cfi::order, rho_, rhoExtra_, rhoM_, findQualityFiles::size, towExcludeEta_, towExcludePhi_, and towExcludePt_.
Referenced by produce().
◆ setupGeometryMap()
Definition at line 236 of file HiPuRhoProducer.cc.
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();
References allgeomid_, caloGeometryToken_, PV3DBase< T, PVType, FrameType >::eta(), geo_, geomtowers_, edm::EventSetup::getData(), CaloGeometry::getPosition(), CaloGeometry::getValidDetIds(), submitPVValidationJobs::gt, DetId::Hcal, HcalDetId::ieta(), ietamax_, ietamin_, HcalDetId::iphi(), LogDebug, nEtaTow_, PV3DBase< T, PVType, FrameType >::phi(), towermap_, and vngeom_.
Referenced by produce().
◆ subtractPedestal()
void HiPuRhoProducer::subtractPedestal |
( |
std::vector< fastjet::PseudoJet > & |
coll | ) |
|
|
virtual |
Definition at line 423 of file HiPuRhoProducer.cc.
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);
References dropZeroTowers_, emean_, esigma_, CaloTower::et(), CaloTower::ieta(), inputs_, and LogDebug.
Referenced by produce().
◆ allgeomid_
std::vector<HcalDetId> HiPuRhoProducer::allgeomid_ |
|
private |
◆ caloGeometryToken_
◆ caloTowerToken_
◆ dropZeroTowers_
const bool HiPuRhoProducer::dropZeroTowers_ |
|
private |
◆ emean_
std::map<int, double> HiPuRhoProducer::emean_ |
|
private |
◆ esigma_
std::map<int, double> HiPuRhoProducer::esigma_ |
|
private |
◆ etaEdgeHi_
std::vector<double> HiPuRhoProducer::etaEdgeHi_ |
|
private |
◆ etaEdgeLow_
std::vector<double> HiPuRhoProducer::etaEdgeLow_ |
|
private |
◆ etaEdges_
std::vector<double> HiPuRhoProducer::etaEdges_ |
|
private |
◆ eTop4_
std::map<int, std::array<double, 4> > HiPuRhoProducer::eTop4_ |
|
private |
◆ fjClusterSeq_
◆ fjInputs_
std::vector<fastjet::PseudoJet> HiPuRhoProducer::fjInputs_ |
|
private |
◆ fjJetDefinition_
◆ fjJets_
std::vector<fastjet::PseudoJet> HiPuRhoProducer::fjJets_ |
|
private |
◆ fjOriginalInputs_
std::vector<fastjet::PseudoJet> HiPuRhoProducer::fjOriginalInputs_ |
|
private |
◆ geo_
◆ geomtowers_
std::map<int, int> HiPuRhoProducer::geomtowers_ |
|
private |
◆ ietamax_
int HiPuRhoProducer::ietamax_ |
|
private |
◆ ietamin_
int HiPuRhoProducer::ietamin_ |
|
private |
◆ initialValue_
const int HiPuRhoProducer::initialValue_ = -99 |
|
private |
◆ inputs_
std::vector<const CaloTower*> HiPuRhoProducer::inputs_ |
|
private |
◆ medianWindowWidth_
const int HiPuRhoProducer::medianWindowWidth_ |
|
private |
◆ minimumTowersFraction_
const double HiPuRhoProducer::minimumTowersFraction_ |
|
private |
◆ nEtaTow_
constexpr static int HiPuRhoProducer::nEtaTow_ = 82 |
|
staticconstexprprivate |
◆ nSigmaPU_
const double HiPuRhoProducer::nSigmaPU_ |
|
private |
◆ nTow_
std::vector<int> HiPuRhoProducer::nTow_ |
|
private |
◆ ntowersWithJets_
std::map<int, int> HiPuRhoProducer::ntowersWithJets_ |
|
private |
◆ postOrphan_
bool HiPuRhoProducer::postOrphan_ |
|
private |
◆ puPtMin_
const double HiPuRhoProducer::puPtMin_ |
|
private |
◆ radiusPU_
const double HiPuRhoProducer::radiusPU_ |
|
private |
◆ rho_
std::vector<double> HiPuRhoProducer::rho_ |
|
private |
◆ rhoExtra_
std::vector<double> HiPuRhoProducer::rhoExtra_ |
|
private |
◆ rhoM_
std::vector<double> HiPuRhoProducer::rhoM_ |
|
private |
◆ rParam_
const double HiPuRhoProducer::rParam_ |
|
private |
◆ setInitialValue_
bool HiPuRhoProducer::setInitialValue_ |
|
private |
◆ towermap_
◆ towExcludeEta_
std::vector<double> HiPuRhoProducer::towExcludeEta_ |
|
private |
◆ towExcludePhi_
std::vector<double> HiPuRhoProducer::towExcludePhi_ |
|
private |
◆ towExcludePt_
std::vector<double> HiPuRhoProducer::towExcludePt_ |
|
private |
◆ towSigmaCut_
const double HiPuRhoProducer::towSigmaCut_ |
|
private |
◆ vmean0_
std::array<float, nEtaTow_> HiPuRhoProducer::vmean0_ |
|
private |
◆ vmean1_
std::array<float, nEtaTow_> HiPuRhoProducer::vmean1_ |
|
private |
◆ vngeom_
std::array<int, nEtaTow_> HiPuRhoProducer::vngeom_ |
|
private |
◆ vntow_
std::array<int, nEtaTow_> HiPuRhoProducer::vntow_ |
|
private |
◆ vrho0_
std::array<float, nEtaTow_> HiPuRhoProducer::vrho0_ |
|
private |
◆ vrho1_
std::array<float, nEtaTow_> HiPuRhoProducer::vrho1_ |
|
private |
◆ vrms0_
std::array<float, nEtaTow_> HiPuRhoProducer::vrms0_ |
|
private |
◆ vrms1_
std::array<float, nEtaTow_> HiPuRhoProducer::vrms1_ |
|
private |
std::vector< fastjet::PseudoJet > fjInputs_
std::vector< double > towExcludePt_
std::vector< double > towExcludeEta_
static const std::string input
std::map< int, double > emean_
constexpr static int nEtaTow_
std::vector< fastjet::PseudoJet > fjJets_
constexpr bool isNotFinite(T x)
constexpr int iphi() const
get the cell iphi
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
std::array< int, nEtaTow_ > vngeom_
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
std::vector< const CaloTower * > inputs_
std::map< int, int > ntowersWithJets_
double pt() const final
transverse momentum
const double towSigmaCut_
ClusterSequencePtr fjClusterSeq_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
virtual void subtractPedestal(std::vector< fastjet::PseudoJet > &coll)
std::vector< double > rhoM_
virtual void calculateOrphanInput(std::vector< fastjet::PseudoJet > &orphanInput)
std::array< float, nEtaTow_ > vrho1_
std::vector< double > etaEdgeLow_
JetDefPtr fjJetDefinition_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double py() const final
y coordinate of momentum vector
constexpr std::array< double, 42 > etaedge
const bool dropZeroTowers_
std::vector< double > towExcludePhi_
std::vector< double > rho_
std::vector< EtaPhiTower > towermap_
virtual void putRho(edm::Event &iEvent, const edm::EventSetup &iSetup)
const double minimumTowersFraction_
constexpr int ieta() const
get the cell ieta
double eta() const final
momentum pseudorapidity
std::array< float, nEtaTow_ > vrms1_
std::array< float, nEtaTow_ > vrho0_
std::vector< fastjet::PseudoJet > fjOriginalInputs_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
std::map< int, std::array< double, 4 > > eTop4_
double et(double vtxZ) const
std::vector< double > rhoExtra_
std::array< int, nEtaTow_ > vntow_
bool getData(T &iHolder) const
virtual void inputTowers()
std::pair< double, double > EtaPhi
std::vector< HcalDetId > allgeomid_
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
const edm::EDGetTokenT< CaloTowerCollection > caloTowerToken_
double phi() const final
momentum azimuthal angle
std::vector< double > etaEdges_
std::array< float, nEtaTow_ > vmean1_
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryToken_
std::array< float, nEtaTow_ > vmean0_
T getParameter(std::string const &) const
double energy() const final
energy
virtual void setupGeometryMap(edm::Event &iEvent, const edm::EventSetup &iSetup)
Abs< T >::type abs(const T &t)
std::array< float, nEtaTow_ > vrms0_
const int medianWindowWidth_
CaloGeometry const * geo_
std::map< int, int > geomtowers_
std::vector< double > etaEdgeHi_
double px() const final
x coordinate of momentum vector
std::map< int, double > esigma_
double pz() const final
z coordinate of momentum vector
Geom::Phi< T > phi() const
virtual void calculatePedestal(std::vector< fastjet::PseudoJet > const &coll)