CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Types | Public Member Functions | Protected Attributes
PileUpSubtractor Class Reference

#include <PileUpSubtractor.h>

Inheritance diagram for PileUpSubtractor:
JetOffsetCorrector MultipleAlgoIterator ParametrizedSubtractor ReflectedIterator

Public Types

typedef std::shared_ptr
< fastjet::GhostedAreaSpec > 
ActiveAreaSpecPtr
 
typedef std::shared_ptr
< fastjet::ClusterSequence > 
ClusterSequencePtr
 
typedef std::shared_ptr
< fastjet::JetDefinition > 
JetDefPtr
 
typedef std::shared_ptr
< fastjet::RangeDefinition > 
RangeDefPtr
 

Public Member Functions

virtual void calculateOrphanInput (std::vector< fastjet::PseudoJet > &orphanInput)
 
virtual void calculatePedestal (std::vector< fastjet::PseudoJet > const &coll)
 
virtual double getCone (double cone, double eta, double phi, double &et, double &pu)
 
virtual double getMeanAtTower (const reco::CandidatePtr &in) const
 
int getN (const reco::CandidatePtr &in) const
 
int getNwithJets (const reco::CandidatePtr &in) const
 
virtual double getPileUpAtTower (const reco::CandidatePtr &in) const
 
virtual double getPileUpEnergy (int ijet) const
 
virtual double getSigmaAtTower (const reco::CandidatePtr &in) const
 
int ieta (const reco::CandidatePtr &in) const
 
int iphi (const reco::CandidatePtr &in) const
 
virtual void offsetCorrectJets ()
 
 PileUpSubtractor (const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
 
virtual void reset (std::vector< edm::Ptr< reco::Candidate > > &input, std::vector< fastjet::PseudoJet > &towers, std::vector< fastjet::PseudoJet > &output)
 
virtual void setDefinition (JetDefPtr const &jetDef)
 
virtual void setupGeometryMap (edm::Event &iEvent, const edm::EventSetup &iSetup)
 
virtual void subtractPedestal (std::vector< fastjet::PseudoJet > &coll)
 
virtual ~PileUpSubtractor ()
 

Protected Attributes

int activeAreaRepeats
 
std::vector< HcalDetIdallgeomid_
 
bool doAreaFastjet_
 
bool doRhoFastjet_
 
std::map< int, double > emean_
 
std::map< int, double > esigma_
 
ActiveAreaSpecPtr fjActiveArea_
 
ClusterSequencePtr fjClusterSeq_
 
std::vector< fastjet::PseudoJet > * fjInputs_
 
JetDefPtr fjJetDefinition_
 
std::vector< fastjet::PseudoJet > * fjJets_
 
std::vector< fastjet::PseudoJet > fjOriginalInputs_
 
CaloGeometry const * geo_
 
std::map< int, int > geomtowers_
 
edm::ESGetToken< CaloGeometry,
CaloGeometryRecord
geoToken_
 
double ghostArea
 
double ghostEtaMax
 
int ietamax_
 
int ietamin_
 
std::vector< edm::Ptr
< reco::Candidate > > * 
inputs_
 
std::vector< double > jetOffset_
 
double jetPtMin_
 
double nSigmaPU_
 
std::map< int, int > ntowersWithJets_
 
double puPtMin_
 
double radiusPU_
 
bool reRunAlgo_
 

Detailed Description

Definition at line 25 of file PileUpSubtractor.h.

Member Typedef Documentation

typedef std::shared_ptr<fastjet::GhostedAreaSpec> PileUpSubtractor::ActiveAreaSpecPtr

Definition at line 28 of file PileUpSubtractor.h.

typedef std::shared_ptr<fastjet::ClusterSequence> PileUpSubtractor::ClusterSequencePtr

Definition at line 27 of file PileUpSubtractor.h.

typedef std::shared_ptr<fastjet::JetDefinition> PileUpSubtractor::JetDefPtr

Definition at line 30 of file PileUpSubtractor.h.

typedef std::shared_ptr<fastjet::RangeDefinition> PileUpSubtractor::RangeDefPtr

Definition at line 29 of file PileUpSubtractor.h.

Constructor & Destructor Documentation

PileUpSubtractor::PileUpSubtractor ( const edm::ParameterSet iConfig,
edm::ConsumesCollector &&  iC 
)

Definition at line 20 of file PileUpSubtractor.cc.

References edm::ParameterSet::getParameter().

20  {
21  geo_ = nullptr;
22  geoToken_ = iC.esConsumes();
23  doAreaFastjet_ = iConfig.getParameter<bool>("doAreaFastjet");
24  doRhoFastjet_ = iConfig.getParameter<bool>("doRhoFastjet");
25  nSigmaPU_ = iConfig.getParameter<double>("nSigmaPU");
26  radiusPU_ = iConfig.getParameter<double>("radiusPU");
27  jetPtMin_ = iConfig.getParameter<double>("jetPtMin");
28  puPtMin_ = iConfig.getParameter<double>("puPtMin");
29  ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
30  activeAreaRepeats = iConfig.getParameter<int>("Active_Area_Repeats");
31  ghostArea = iConfig.getParameter<double>("GhostArea");
32 
33  if (doAreaFastjet_ || doRhoFastjet_) {
34  fjActiveArea_ = std::make_shared<fastjet::ActiveAreaSpec>(ghostEtaMax, activeAreaRepeats, ghostArea);
35  if ((ghostEtaMax < 0) || (activeAreaRepeats < 0) || (ghostArea < 0))
36  throw cms::Exception("doAreaFastjet or doRhoFastjet")
37  << "Parameters ghostEtaMax, activeAreaRepeats or ghostArea for doAreaFastjet/doRhoFastjet are not defined."
38  << std::endl;
39  }
40 }
ActiveAreaSpecPtr fjActiveArea_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
CaloGeometry const * geo_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geoToken_
virtual PileUpSubtractor::~PileUpSubtractor ( )
inlinevirtual

Definition at line 33 of file PileUpSubtractor.h.

33 { ; }

Member Function Documentation

void PileUpSubtractor::calculateOrphanInput ( std::vector< fastjet::PseudoJet > &  orphanInput)
virtual

Reimplemented in ParametrizedSubtractor, and MultipleAlgoIterator.

Definition at line 183 of file PileUpSubtractor.cc.

References reco::deltaR(), runTauDisplay::dr, spr::find(), and LogDebug.

183  {
184  LogDebug("PileUpSubtractor") << "The subtractor calculating orphan input...\n";
185 
186  (*fjInputs_) = fjOriginalInputs_;
187 
188  vector<int> jettowers; // vector of towers indexed by "user_index"
189  vector<pair<int, int> > excludedTowers; // vector of excluded ieta, iphi values
190 
191  vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin(), fjJetsEnd = fjJets_->end();
192  for (; pseudojetTMP != fjJetsEnd; ++pseudojetTMP) {
193  if (pseudojetTMP->perp() < puPtMin_)
194  continue;
195 
196  // find towers within radiusPU_ of this jet
197  for (vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++) {
198  double dr = reco::deltaR(geo_->getPosition((DetId)(*im)), (*pseudojetTMP));
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()));
204  }
205  }
206  vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(), fjInputsEnd = fjInputs_->end();
207 
208  for (; it != fjInputsEnd; ++it) {
209  int index = it->user_index();
210  int ie = ieta((*inputs_)[index]);
211  int ip = iphi((*inputs_)[index]);
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);
216  } //dr < radiusPU_
217  } // initial input collection
218  } // pseudojets
219 
220  //
221  // Create a new collections from the towers not included in jets
222  //
223  for (vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(), fjInputsEnd = fjInputs_->end();
224  it != fjInputsEnd;
225  ++it) {
226  int index = it->user_index();
227  vector<int>::const_iterator itjet = find(jettowers.begin(), jettowers.end(), index);
228  if (itjet == jettowers.end()) {
229  const reco::CandidatePtr& originalTower = (*inputs_)[index];
230  fastjet::PseudoJet orphan(originalTower->px(), originalTower->py(), originalTower->pz(), originalTower->energy());
231  orphan.set_user_index(index);
232 
233  orphanInput.push_back(orphan);
234  }
235  }
236 }
std::vector< fastjet::PseudoJet > * fjJets_
std::vector< fastjet::PseudoJet > fjOriginalInputs_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
int ieta(const reco::CandidatePtr &in) const
std::vector< fastjet::PseudoJet > * fjInputs_
std::map< int, int > ntowersWithJets_
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
Definition: DetId.h:17
int iphi(const reco::CandidatePtr &in) const
CaloGeometry const * geo_
std::vector< HcalDetId > allgeomid_
std::vector< edm::Ptr< reco::Candidate > > * inputs_
#define LogDebug(id)
void PileUpSubtractor::calculatePedestal ( std::vector< fastjet::PseudoJet > const &  coll)
virtual

Reimplemented in ParametrizedSubtractor, ReflectedIterator, and MultipleAlgoIterator.

Definition at line 95 of file PileUpSubtractor.cc.

References submitPVValidationJobs::gt, mps_fire::i, LogDebug, nt, edm::second(), and mathSSE::sqrt().

95  {
96  LogDebug("PileUpSubtractor") << "The subtractor calculating pedestals...\n";
97  map<int, double> emean2;
98  map<int, int> ntowers;
99 
100  int ietaold = -10000;
101  int ieta0 = -100;
102 
103  // Initial values for emean_, emean2, esigma_, ntowers
104 
105  for (int i = ietamin_; i < ietamax_ + 1; i++) {
106  emean_[i] = 0.;
107  emean2[i] = 0.;
108  esigma_[i] = 0.;
109  ntowers[i] = 0;
110  }
111 
112  for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin(), fjInputsEnd = coll.end();
113  input_object != fjInputsEnd;
114  ++input_object) {
115  const reco::CandidatePtr& originalTower = (*inputs_)[input_object->user_index()];
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;
121  ntowers[ieta0] = 1;
122  ietaold = ieta0;
123  } else {
124  emean_[ieta0] = emean_[ieta0] + Original_Et;
125  emean2[ieta0] = emean2[ieta0] + Original_Et * Original_Et;
126  ntowers[ieta0]++;
127  }
128  }
129 
130  for (map<int, int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++) {
131  int it = (*gt).first;
132 
133  double e1 = (*(emean_.find(it))).second;
134  double e2 = (*emean2.find(it)).second;
135  int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
136 
137  LogDebug("PileUpSubtractor") << " ieta : " << it << " number of towers : " << nt << " e1 : " << e1 << " e2 : " << e2
138  << "\n";
139  if (nt > 0) {
140  emean_[it] = e1 / nt;
141  double eee = e2 / nt - e1 * e1 / (nt * nt);
142  if (eee < 0.)
143  eee = 0.;
144  esigma_[it] = nSigmaPU_ * sqrt(eee);
145  } else {
146  emean_[it] = 0.;
147  esigma_[it] = 0.;
148  }
149  LogDebug("PileUpSubtractor") << " ieta : " << it << " Pedestals : " << emean_[it] << " " << esigma_[it] << "\n";
150  }
151 }
std::map< int, double > esigma_
int ieta(const reco::CandidatePtr &in) const
std::map< int, int > geomtowers_
U second(std::pair< T, U > const &p)
std::map< int, double > emean_
std::map< int, int > ntowersWithJets_
T sqrt(T t)
Definition: SSEVec.h:19
int nt
Definition: AMPTWrapper.h:42
#define LogDebug(id)
double PileUpSubtractor::getCone ( double  cone,
double  eta,
double  phi,
double &  et,
double &  pu 
)
virtual

Definition at line 280 of file PileUpSubtractor.cc.

References reco::deltaR(), runTauDisplay::dr, PVValHelper::eta, PV3DBase< T, PVType, FrameType >::eta(), PV3DBase< T, PVType, FrameType >::phi(), and point.

280  {
281  pu = 0;
282 
283  for (vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++) {
284  if (im->depth() != 1)
285  continue;
286  const GlobalPoint& point = geo_->getPosition((DetId)(*im));
287  double dr = reco::deltaR(point.eta(), point.phi(), eta, phi);
288  if (dr < cone) {
289  pu += (*emean_.find(im->ieta())).second + (*esigma_.find(im->ieta())).second;
290  }
291  }
292 
293  return pu;
294 }
std::map< int, double > esigma_
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
std::map< int, double > emean_
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
Definition: DetId.h:17
T eta() const
Definition: PV3DBase.h:73
CaloGeometry const * geo_
std::vector< HcalDetId > allgeomid_
*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
Definition: invegas.h:5
double PileUpSubtractor::getMeanAtTower ( const reco::CandidatePtr in) const
virtual

Reimplemented in ParametrizedSubtractor.

Definition at line 296 of file PileUpSubtractor.cc.

References edm::second().

296  {
297  int it = ieta(in);
298  return (*emean_.find(it)).second;
299 }
int ieta(const reco::CandidatePtr &in) const
U second(std::pair< T, U > const &p)
std::map< int, double > emean_
int PileUpSubtractor::getN ( const reco::CandidatePtr in) const

Definition at line 311 of file PileUpSubtractor.cc.

References dqmiodumpmetadata::n.

311  {
312  int it = ieta(in);
313 
314  int n = (*(geomtowers_.find(it))).second - (*(ntowersWithJets_.find(it))).second;
315  return n;
316 }
int ieta(const reco::CandidatePtr &in) const
std::map< int, int > geomtowers_
std::map< int, int > ntowersWithJets_
int PileUpSubtractor::getNwithJets ( const reco::CandidatePtr in) const

Definition at line 318 of file PileUpSubtractor.cc.

References dqmiodumpmetadata::n.

318  {
319  int it = ieta(in);
320  int n = (*(ntowersWithJets_.find(it))).second;
321  return n;
322 }
int ieta(const reco::CandidatePtr &in) const
std::map< int, int > ntowersWithJets_
double PileUpSubtractor::getPileUpAtTower ( const reco::CandidatePtr in) const
virtual

Reimplemented in ParametrizedSubtractor.

Definition at line 306 of file PileUpSubtractor.cc.

References edm::second().

306  {
307  int it = ieta(in);
308  return (*emean_.find(it)).second + (*esigma_.find(it)).second;
309 }
std::map< int, double > esigma_
int ieta(const reco::CandidatePtr &in) const
U second(std::pair< T, U > const &p)
std::map< int, double > emean_
virtual double PileUpSubtractor::getPileUpEnergy ( int  ijet) const
inlinevirtual

Definition at line 47 of file PileUpSubtractor.h.

References jetOffset_.

47 { return jetOffset_[ijet]; }
std::vector< double > jetOffset_
double PileUpSubtractor::getSigmaAtTower ( const reco::CandidatePtr in) const
virtual

Reimplemented in ParametrizedSubtractor.

Definition at line 301 of file PileUpSubtractor.cc.

References edm::second().

301  {
302  int it = ieta(in);
303  return (*esigma_.find(it)).second;
304 }
std::map< int, double > esigma_
int ieta(const reco::CandidatePtr &in) const
U second(std::pair< T, U > const &p)
int PileUpSubtractor::ieta ( const reco::CandidatePtr in) const

Definition at line 324 of file PileUpSubtractor.cc.

References Exception, edm::Ptr< T >::get(), CaloTower::id(), and CaloTowerDetId::ieta().

Referenced by MultipleAlgoIterator::calculateOrphanInput(), MultipleAlgoIterator::calculatePedestal(), ParametrizedSubtractor::getMeanAtTower(), ParametrizedSubtractor::getPileUpAtTower(), ParametrizedSubtractor::getSigmaAtTower(), MultipleAlgoIterator::offsetCorrectJets(), ParametrizedSubtractor::offsetCorrectJets(), MultipleAlgoIterator::subtractPedestal(), and ParametrizedSubtractor::subtractPedestal().

324  {
325  int it = 0;
326  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
327  if (ctc) {
328  it = ctc->id().ieta();
329  } else {
330  throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
331  }
332  return it;
333 }
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:139
CaloTowerDetId id() const
Definition: CaloTower.h:123
int ieta() const
get the tower ieta
int PileUpSubtractor::iphi ( const reco::CandidatePtr in) const

Definition at line 335 of file PileUpSubtractor.cc.

References Exception, edm::Ptr< T >::get(), CaloTower::id(), and CaloTowerDetId::iphi().

Referenced by MultipleAlgoIterator::calculateOrphanInput(), MultipleAlgoIterator::subtractPedestal(), and ParametrizedSubtractor::subtractPedestal().

335  {
336  int it = 0;
337  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
338  if (ctc) {
339  it = ctc->id().iphi();
340  } else {
341  throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
342  }
343  return it;
344 }
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:139
int iphi() const
get the tower iphi
CaloTowerDetId id() const
Definition: CaloTower.h:123
void PileUpSubtractor::offsetCorrectJets ( )
virtual

Reimplemented in ParametrizedSubtractor, ReflectedIterator, and MultipleAlgoIterator.

Definition at line 238 of file PileUpSubtractor.cc.

References LogDebug, dt_dqm_sourceclient_common_cff::reco, edm::second(), and HLT_FULL_cff::towers.

238  {
239  LogDebug("PileUpSubtractor") << "The subtractor correcting jets...\n";
240  jetOffset_.clear();
241  using namespace reco;
242 
243  //
244  // Reestimate energy of jet (energy of jet with initial map)
245  //
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;
251 
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) {
255  const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
256  int it = ieta(originalTower);
257  double Original_Et = originalTower->et();
258  double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second;
259  if (etnew < 0.)
260  etnew = 0;
261  newjetet = newjetet + etnew;
262  jetOffset_[ijet] += Original_Et - etnew;
263  }
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]
269  << '\n';
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);
277  }
278 }
std::vector< double > jetOffset_
std::map< int, double > esigma_
std::vector< fastjet::PseudoJet > * fjJets_
int ieta(const reco::CandidatePtr &in) const
U second(std::pair< T, U > const &p)
std::map< int, double > emean_
#define LogDebug(id)
void PileUpSubtractor::reset ( std::vector< edm::Ptr< reco::Candidate > > &  input,
std::vector< fastjet::PseudoJet > &  towers,
std::vector< fastjet::PseudoJet > &  output 
)
virtual

Definition at line 42 of file PileUpSubtractor.cc.

References mps_fire::i, input, convertSQLitetoXML_cfg::output, and HLT_FULL_cff::towers.

44  {
45  inputs_ = &input;
46  fjInputs_ = &towers;
47  fjJets_ = &output;
48  fjOriginalInputs_ = (*fjInputs_);
49  for (unsigned int i = 0; i < fjInputs_->size(); ++i) {
50  fjOriginalInputs_[i].set_user_index((*fjInputs_)[i].user_index());
51  }
52 }
std::vector< fastjet::PseudoJet > * fjJets_
std::vector< fastjet::PseudoJet > fjOriginalInputs_
static std::string const input
Definition: EdmProvDump.cc:47
std::vector< fastjet::PseudoJet > * fjInputs_
std::vector< edm::Ptr< reco::Candidate > > * inputs_
void PileUpSubtractor::setDefinition ( JetDefPtr const &  jetDef)
virtual

Definition at line 54 of file PileUpSubtractor.cc.

54  {
55  fjJetDefinition_ = std::make_shared<fastjet::JetDefinition>(*jetDef);
56 }
JetDefPtr fjJetDefinition_
void PileUpSubtractor::setupGeometryMap ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Reimplemented in ParametrizedSubtractor.

Definition at line 58 of file PileUpSubtractor.cc.

References edm::EventSetup::getData(), DetId::Hcal, mps_fire::i, HcalDetId::ieta(), and LogDebug.

Referenced by ParametrizedSubtractor::setupGeometryMap().

58  {
59  LogDebug("PileUpSubtractor") << "The subtractor setting up geometry...\n";
60 
61  if (geo_ == nullptr) {
62  geo_ = &iSetup.getData(geoToken_);
63  std::vector<DetId> alldid = geo_->getValidDetIds();
64 
65  int ietaold = -10000;
66  ietamax_ = -10000;
67  ietamin_ = 10000;
68  for (std::vector<DetId>::const_iterator did = alldid.begin(); did != alldid.end(); did++) {
69  if ((*did).det() == DetId::Hcal) {
70  HcalDetId hid = HcalDetId(*did);
71  allgeomid_.push_back(*did);
72 
73  if (hid.ieta() != ietaold) {
74  ietaold = hid.ieta();
75  geomtowers_[hid.ieta()] = 1;
76  if (hid.ieta() > ietamax_)
77  ietamax_ = hid.ieta();
78  if (hid.ieta() < ietamin_)
79  ietamin_ = hid.ieta();
80  } else {
81  geomtowers_[hid.ieta()]++;
82  }
83  }
84  }
85  }
86 
87  for (int i = ietamin_; i < ietamax_ + 1; i++) {
88  ntowersWithJets_[i] = 0;
89  }
90 }
std::map< int, int > geomtowers_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
std::map< int, int > ntowersWithJets_
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:75
CaloGeometry const * geo_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geoToken_
std::vector< HcalDetId > allgeomid_
#define LogDebug(id)
void PileUpSubtractor::subtractPedestal ( std::vector< fastjet::PseudoJet > &  coll)
virtual

Reimplemented in ParametrizedSubtractor, ReflectedIterator, and MultipleAlgoIterator.

Definition at line 156 of file PileUpSubtractor.cc.

References LogDebug.

156  {
157  LogDebug("PileUpSubtractor") << "The subtractor subtracting pedestals...\n";
158 
159  int it = -100;
160  for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin(), fjInputsEnd = coll.end();
161  input_object != fjInputsEnd;
162  ++input_object) {
163  reco::CandidatePtr const& itow = (*inputs_)[input_object->user_index()];
164 
165  it = ieta(itow);
166 
167  double etnew = itow->et() - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
168  float mScale = etnew / input_object->Et();
169  if (etnew < 0.)
170  mScale = 0.;
171 
172  math::XYZTLorentzVectorD towP4(input_object->px() * mScale,
173  input_object->py() * mScale,
174  input_object->pz() * mScale,
175  input_object->e() * mScale);
176 
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);
180  }
181 }
std::map< int, double > esigma_
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
int ieta(const reco::CandidatePtr &in) const
std::map< int, double > emean_
#define LogDebug(id)

Member Data Documentation

int PileUpSubtractor::activeAreaRepeats
protected

Definition at line 72 of file PileUpSubtractor.h.

std::vector<HcalDetId> PileUpSubtractor::allgeomid_
protected

Definition at line 82 of file PileUpSubtractor.h.

Referenced by MultipleAlgoIterator::calculateOrphanInput().

bool PileUpSubtractor::doAreaFastjet_
protected
bool PileUpSubtractor::doRhoFastjet_
protected
std::map<int, double> PileUpSubtractor::emean_
protected
std::map<int, double> PileUpSubtractor::esigma_
protected
ActiveAreaSpecPtr PileUpSubtractor::fjActiveArea_
protected
ClusterSequencePtr PileUpSubtractor::fjClusterSeq_
protected
std::vector<fastjet::PseudoJet>* PileUpSubtractor::fjInputs_
protected
JetDefPtr PileUpSubtractor::fjJetDefinition_
protected

Definition at line 57 of file PileUpSubtractor.h.

Referenced by MultipleAlgoIterator::offsetCorrectJets().

std::vector<fastjet::PseudoJet>* PileUpSubtractor::fjJets_
protected
std::vector<fastjet::PseudoJet> PileUpSubtractor::fjOriginalInputs_
protected
CaloGeometry const* PileUpSubtractor::geo_
protected
std::map<int, int> PileUpSubtractor::geomtowers_
protected
edm::ESGetToken<CaloGeometry, CaloGeometryRecord> PileUpSubtractor::geoToken_
protected

Definition at line 78 of file PileUpSubtractor.h.

double PileUpSubtractor::ghostArea
protected

Definition at line 73 of file PileUpSubtractor.h.

double PileUpSubtractor::ghostEtaMax
protected

Definition at line 71 of file PileUpSubtractor.h.

int PileUpSubtractor::ietamax_
protected
int PileUpSubtractor::ietamin_
protected
std::vector<edm::Ptr<reco::Candidate> >* PileUpSubtractor::inputs_
protected

Definition at line 59 of file PileUpSubtractor.h.

Referenced by MultipleAlgoIterator::calculateOrphanInput().

std::vector<double> PileUpSubtractor::jetOffset_
protected
double PileUpSubtractor::jetPtMin_
protected
double PileUpSubtractor::nSigmaPU_
protected
std::map<int, int> PileUpSubtractor::ntowersWithJets_
protected
double PileUpSubtractor::puPtMin_
protected

Definition at line 69 of file PileUpSubtractor.h.

Referenced by MultipleAlgoIterator::calculateOrphanInput().

double PileUpSubtractor::radiusPU_
protected

Definition at line 76 of file PileUpSubtractor.h.

Referenced by MultipleAlgoIterator::calculateOrphanInput().

bool PileUpSubtractor::reRunAlgo_
protected

Definition at line 65 of file PileUpSubtractor.h.