22 geoToken_ = iC.esConsumes();
23 doAreaFastjet_ = iConfig.
getParameter<
bool>(
"doAreaFastjet");
24 doRhoFastjet_ = iConfig.
getParameter<
bool>(
"doRhoFastjet");
29 ghostEtaMax = iConfig.
getParameter<
double>(
"Ghost_EtaMax");
30 activeAreaRepeats = iConfig.
getParameter<
int>(
"Active_Area_Repeats");
33 if (doAreaFastjet_ || doRhoFastjet_) {
34 fjActiveArea_ = std::make_shared<fastjet::ActiveAreaSpec>(ghostEtaMax, activeAreaRepeats, ghostArea);
35 if ((ghostEtaMax < 0) || (activeAreaRepeats < 0) || (ghostArea < 0))
37 <<
"Parameters ghostEtaMax, activeAreaRepeats or ghostArea for doAreaFastjet/doRhoFastjet are not defined." 43 std::vector<fastjet::PseudoJet>&
towers,
44 std::vector<fastjet::PseudoJet>&
output) {
48 fjOriginalInputs_ = (*fjInputs_);
49 for (
unsigned int i = 0;
i < fjInputs_->size(); ++
i) {
50 fjOriginalInputs_[
i].set_user_index((*fjInputs_)[
i].user_index());
55 fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(*jetDef);
59 LogDebug(
"PileUpSubtractor") <<
"The subtractor setting up geometry...\n";
61 if (geo_ ==
nullptr) {
62 geo_ = &iSetup.
getData(geoToken_);
63 std::vector<DetId> alldid = geo_->getValidDetIds();
68 for (std::vector<DetId>::const_iterator did = alldid.begin(); did != alldid.end(); did++) {
71 allgeomid_.push_back(*did);
73 if (hid.
ieta() != ietaold) {
75 geomtowers_[hid.
ieta()] = 1;
76 if (hid.
ieta() > ietamax_)
77 ietamax_ = hid.
ieta();
78 if (hid.
ieta() < ietamin_)
79 ietamin_ = hid.
ieta();
81 geomtowers_[hid.
ieta()]++;
87 for (
int i = ietamin_;
i < ietamax_ + 1;
i++) {
88 ntowersWithJets_[
i] = 0;
96 LogDebug(
"PileUpSubtractor") <<
"The subtractor calculating pedestals...\n";
97 map<int, double> emean2;
98 map<int, int> ntowers;
100 int ietaold = -10000;
105 for (
int i = ietamin_;
i < ietamax_ + 1;
i++) {
112 for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin(), fjInputsEnd = coll.end();
113 input_object != fjInputsEnd;
116 ieta0 =
ieta(originalTower);
117 double Original_Et = originalTower->et();
118 if (ieta0 - ietaold != 0) {
119 emean_[ieta0] = emean_[ieta0] + Original_Et;
120 emean2[ieta0] = emean2[ieta0] + Original_Et * Original_Et;
124 emean_[ieta0] = emean_[ieta0] + Original_Et;
125 emean2[ieta0] = emean2[ieta0] + Original_Et * Original_Et;
130 for (map<int, int>::const_iterator
gt = geomtowers_.begin();
gt != geomtowers_.end();
gt++) {
131 int it = (*gt).first;
133 double e1 = (*(emean_.find(
it))).second;
134 double e2 = (*emean2.find(
it)).
second;
135 int nt = (*gt).second - (*(ntowersWithJets_.find(
it))).second;
137 LogDebug(
"PileUpSubtractor") <<
" ieta : " <<
it <<
" number of towers : " <<
nt <<
" e1 : " <<
e1 <<
" e2 : " << e2
144 esigma_[
it] = nSigmaPU_ *
sqrt(eee);
149 LogDebug(
"PileUpSubtractor") <<
" ieta : " <<
it <<
" Pedestals : " << emean_[
it] <<
" " << esigma_[
it] <<
"\n";
157 LogDebug(
"PileUpSubtractor") <<
"The subtractor subtracting pedestals...\n";
160 for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin(), fjInputsEnd = coll.end();
161 input_object != fjInputsEnd;
167 double etnew = itow->et() - (*(emean_.find(
it))).second - (*(esigma_.find(
it))).second;
168 float mScale = etnew / input_object->Et();
173 input_object->py() * mScale,
174 input_object->pz() * mScale,
175 input_object->e() * mScale);
177 int index = input_object->user_index();
178 input_object->reset_momentum(towP4.px(), towP4.py(), towP4.pz(), towP4.energy());
179 input_object->set_user_index(
index);
184 LogDebug(
"PileUpSubtractor") <<
"The subtractor calculating orphan input...\n";
186 (*fjInputs_) = fjOriginalInputs_;
188 vector<int> jettowers;
189 vector<pair<int, int> > excludedTowers;
191 vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin(), fjJetsEnd = fjJets_->end();
192 for (; pseudojetTMP != fjJetsEnd; ++pseudojetTMP) {
193 if (pseudojetTMP->perp() < puPtMin_)
197 for (vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++) {
199 vector<pair<int, int> >::const_iterator exclude =
200 find(excludedTowers.begin(), excludedTowers.end(), pair<int, int>(im->ieta(), im->iphi()));
201 if (dr < radiusPU_ && exclude == excludedTowers.end()) {
202 ntowersWithJets_[(*im).ieta()]++;
203 excludedTowers.push_back(pair<int, int>(im->ieta(), im->iphi()));
206 vector<fastjet::PseudoJet>::const_iterator
it = fjInputs_->begin(), fjInputsEnd = fjInputs_->end();
208 for (;
it != fjInputsEnd; ++
it) {
212 vector<pair<int, int> >::const_iterator exclude =
213 find(excludedTowers.begin(), excludedTowers.end(), pair<int, int>(ie, ip));
214 if (exclude != excludedTowers.end()) {
215 jettowers.push_back(
index);
223 for (vector<fastjet::PseudoJet>::const_iterator
it = fjInputs_->begin(), fjInputsEnd = fjInputs_->end();
227 vector<int>::const_iterator itjet =
find(jettowers.begin(), jettowers.end(),
index);
228 if (itjet == jettowers.end()) {
230 fastjet::PseudoJet orphan(originalTower->px(), originalTower->py(), originalTower->pz(), originalTower->energy());
231 orphan.set_user_index(
index);
233 orphanInput.push_back(orphan);
239 LogDebug(
"PileUpSubtractor") <<
"The subtractor correcting jets...\n";
241 using namespace reco;
246 jetOffset_.reserve(fjJets_->size());
247 vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin(), jetsEnd = fjJets_->end();
248 for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
249 int ijet = pseudojetTMP - fjJets_->begin();
250 jetOffset_[ijet] = 0;
252 std::vector<fastjet::PseudoJet>
towers = fastjet::sorted_by_pt(pseudojetTMP->constituents());
253 double newjetet = 0.;
254 for (vector<fastjet::PseudoJet>::const_iterator ito =
towers.begin(), towEnd =
towers.end(); ito != towEnd; ++ito) {
256 int it =
ieta(originalTower);
257 double Original_Et = originalTower->et();
258 double etnew = Original_Et - (*emean_.find(
it)).
second - (*esigma_.find(
it)).
second;
261 newjetet = newjetet + etnew;
262 jetOffset_[ijet] += Original_Et - etnew;
264 double mScale = newjetet / pseudojetTMP->Et();
265 LogDebug(
"PileUpSubtractor") <<
"pseudojetTMP->Et() : " << pseudojetTMP->Et() <<
'\n';
266 LogDebug(
"PileUpSubtractor") <<
"newjetet : " << newjetet <<
'\n';
267 LogDebug(
"PileUpSubtractor") <<
"jetOffset_[ijet] : " << jetOffset_[ijet] <<
'\n';
268 LogDebug(
"PileUpSubtractor") <<
"pseudojetTMP->Et() - jetOffset_[ijet] : " << pseudojetTMP->Et() - jetOffset_[ijet]
270 LogDebug(
"PileUpSubtractor") <<
"Scale is : " << mScale <<
'\n';
271 int cshist = pseudojetTMP->cluster_hist_index();
272 pseudojetTMP->reset_momentum(pseudojetTMP->px() * mScale,
273 pseudojetTMP->py() * mScale,
274 pseudojetTMP->pz() * mScale,
275 pseudojetTMP->e() * mScale);
276 pseudojetTMP->set_cluster_hist_index(cshist);
283 for (vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++) {
284 if (im->depth() != 1)
289 pu += (*emean_.find(im->ieta())).second + (*esigma_.find(im->ieta())).second;
314 int n = (*(geomtowers_.find(
it))).second - (*(ntowersWithJets_.find(
it))).second;
320 int n = (*(ntowersWithJets_.find(
it))).second;
330 throw cms::Exception(
"Invalid Constituent") <<
"CaloJet constituent is not of CaloTower type";
341 throw cms::Exception(
"Invalid Constituent") <<
"CaloJet constituent is not of CaloTower type";
virtual void calculateOrphanInput(std::vector< fastjet::PseudoJet > &orphanInput)
T getParameter(std::string const &) const
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
virtual double getPileUpAtTower(const reco::CandidatePtr &in) const
int getN(const reco::CandidatePtr &in) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
static std::string const input
U second(std::pair< T, U > const &p)
virtual void offsetCorrectJets()
virtual double getMeanAtTower(const reco::CandidatePtr &in) const
constexpr int ieta() const
get the cell ieta
int getNwithJets(const reco::CandidatePtr &in) const
int iphi(const reco::CandidatePtr &in) const
virtual void subtractPedestal(std::vector< fastjet::PseudoJet > &coll)
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
#define EDM_REGISTER_PLUGINFACTORY(_factory_, _category_)
int iphi() const
get the tower iphi
virtual void setDefinition(JetDefPtr const &jetDef)
int ieta() const
get the tower ieta
int ieta(const reco::CandidatePtr &in) const
PileUpSubtractor(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
virtual void calculatePedestal(std::vector< fastjet::PseudoJet > const &coll)
virtual void setupGeometryMap(edm::Event &iEvent, const edm::EventSetup &iSetup)
virtual double getSigmaAtTower(const reco::CandidatePtr &in) const
virtual double getCone(double cone, double eta, double phi, double &et, double &pu)
virtual void reset(std::vector< edm::Ptr< reco::Candidate > > &input, std::vector< fastjet::PseudoJet > &towers, std::vector< fastjet::PseudoJet > &output)
std::shared_ptr< fastjet::JetDefinition > JetDefPtr
CaloTowerDetId id() const
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point