CMS 3D CMS Logo

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

#include <ReflectedIterator.h>

Inheritance diagram for ReflectedIterator:
PileUpSubtractor

Public Member Functions

void calculatePedestal (std::vector< fastjet::PseudoJet > const &coll) override
 
double getEt (const reco::CandidatePtr &in) const
 
double getEta (const reco::CandidatePtr &in) const
 
void offsetCorrectJets () override
 
 ReflectedIterator (const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
 
void rescaleRMS (double s)
 
void subtractPedestal (std::vector< fastjet::PseudoJet > &coll) override
 
 ~ReflectedIterator () override
 
- Public Member Functions inherited from PileUpSubtractor
virtual void calculateOrphanInput (std::vector< fastjet::PseudoJet > &orphanInput)
 
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
 
 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 ~PileUpSubtractor ()
 

Public Attributes

bool dropZeroTowers_
 
bool sumRecHits_
 

Additional Inherited Members

- Public Types inherited from PileUpSubtractor
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
 
- Protected Attributes inherited from PileUpSubtractor
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 6 of file ReflectedIterator.h.

Constructor & Destructor Documentation

ReflectedIterator::ReflectedIterator ( const edm::ParameterSet iConfig,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 8 of file ReflectedIterator.h.

9  : PileUpSubtractor(iConfig, std::move(iC)),
10  sumRecHits_(iConfig.getParameter<bool>("sumRecHits")),
11  dropZeroTowers_(iConfig.getUntrackedParameter<bool>("dropZeroTowers", true)) {
12  ;
13  }
T getUntrackedParameter(std::string const &, T const &) const
def move
Definition: eostools.py:511
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
PileUpSubtractor(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
ReflectedIterator::~ReflectedIterator ( )
inlineoverride

Definition at line 23 of file ReflectedIterator.h.

23 { ; }

Member Function Documentation

void ReflectedIterator::calculatePedestal ( std::vector< fastjet::PseudoJet > const &  coll)
overridevirtual

Reimplemented from PileUpSubtractor.

Definition at line 99 of file ReflectedIterator.cc.

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

99  {
100  LogDebug("PileUpSubtractor") << "The subtractor calculating pedestals...\n";
101 
102  map<int, double> emean2;
103  map<int, int> ntowers;
104 
105  int ietaold = -10000;
106  int ieta0 = -100;
107 
108  // Initial values for emean_, emean2, esigma_, ntowers
109 
110  for (int i = ietamin_; i < ietamax_ + 1; i++) {
111  emean_[i] = 0.;
112  emean2[i] = 0.;
113  esigma_[i] = 0.;
114  ntowers[i] = 0;
115  }
116 
117  for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin(), fjInputsEnd = coll.end();
118  input_object != fjInputsEnd;
119  ++input_object) {
120  const reco::CandidatePtr& originalTower = (*inputs_)[input_object->user_index()];
121  ieta0 = ieta(originalTower);
122  double Original_Et = originalTower->et();
123  if (sumRecHits_) {
124  Original_Et = getEt(originalTower);
125  }
126 
127  if (ieta0 - ietaold != 0) {
128  emean_[ieta0] = emean_[ieta0] + Original_Et;
129  emean2[ieta0] = emean2[ieta0] + Original_Et * Original_Et;
130  ntowers[ieta0] = 1;
131  ietaold = ieta0;
132  } else {
133  emean_[ieta0] = emean_[ieta0] + Original_Et;
134  emean2[ieta0] = emean2[ieta0] + Original_Et * Original_Et;
135  ntowers[ieta0]++;
136  }
137  }
138 
139  for (map<int, int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++) {
140  int it = (*gt).first;
141 
142  double e1 = (*(emean_.find(it))).second;
143  double e2 = (*emean2.find(it)).second;
144  int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
145 
146  LogDebug("PileUpSubtractor") << " ieta : " << it << " number of towers : " << nt << " e1 : " << e1 << " e2 : " << e2
147  << "\n";
148 
149  if (nt > 0) {
150  emean_[it] = e1 / nt;
151  double eee = e2 / nt - e1 * e1 / (nt * nt);
152  if (eee < 0.)
153  eee = 0.;
154  esigma_[it] = nSigmaPU_ * sqrt(eee);
155  } else {
156  emean_[it] = 0.;
157  esigma_[it] = 0.;
158  }
159  LogDebug("PileUpSubtractor") << " ieta : " << it << " Pedestals : " << emean_[it] << " " << esigma_[it] << "\n";
160  }
161 }
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
double getEt(const reco::CandidatePtr &in) const
#define LogDebug(id)
double ReflectedIterator::getEt ( const reco::CandidatePtr in) const

Definition at line 163 of file ReflectedIterator.cc.

References CaloTower::emEnergy(), relval_parameters_module::energy, edm::Ptr< T >::get(), CaloTower::hadEnergy(), CaloTower::id(), funct::sin(), and PV3DBase< T, PVType, FrameType >::theta().

163  {
164  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
165  const GlobalPoint& pos = geo_->getPosition(ctc->id());
166  double energy = ctc->emEnergy() + ctc->hadEnergy();
167  double et = energy * sin(pos.theta());
168  return et;
169 }
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:139
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
double emEnergy() const
Definition: CaloTower.h:130
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
double hadEnergy() const
Definition: CaloTower.h:131
CaloTowerDetId id() const
Definition: CaloTower.h:123
CaloGeometry const * geo_
double ReflectedIterator::getEta ( const reco::CandidatePtr in) const

Definition at line 171 of file ReflectedIterator.cc.

References PVValHelper::eta, PV3DBase< T, PVType, FrameType >::eta(), edm::Ptr< T >::get(), and CaloTower::id().

171  {
172  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
173  const GlobalPoint& pos = geo_->getPosition(ctc->id());
174  double eta = pos.eta();
175  return eta;
176 }
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:139
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
CaloTowerDetId id() const
Definition: CaloTower.h:123
T eta() const
Definition: PV3DBase.h:73
CaloGeometry const * geo_
void ReflectedIterator::offsetCorrectJets ( )
overridevirtual

Reimplemented from PileUpSubtractor.

Definition at line 14 of file ReflectedIterator.cc.

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

14  {
15  LogDebug("PileUpSubtractor") << "The subtractor correcting jets...\n";
16  jetOffset_.clear();
17 
18  using namespace reco;
19 
20  (*fjInputs_) = fjOriginalInputs_;
23  const fastjet::JetDefinition& def = fjClusterSeq_->jet_def();
24  if (!doAreaFastjet_ && !doRhoFastjet_) {
25  fastjet::ClusterSequence newseq(*fjInputs_, def);
26  (*fjClusterSeq_) = newseq;
27  } else {
28  fastjet::ClusterSequenceArea newseq(*fjInputs_, def, *fjActiveArea_);
29  (*fjClusterSeq_) = newseq;
30  }
31 
32  (*fjJets_) = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
33 
34  jetOffset_.reserve(fjJets_->size());
35 
36  vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin(), jetsEnd = fjJets_->end();
37  for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
38  int ijet = pseudojetTMP - fjJets_->begin();
39  jetOffset_[ijet] = 0;
40 
41  std::vector<fastjet::PseudoJet> towers = sorted_by_pt(fjClusterSeq_->constituents(*pseudojetTMP));
42 
43  double newjetet = 0.;
44  for (vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(), towEnd = towers.end(); ito != towEnd; ++ito) {
45  const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
46  int it = ieta(originalTower);
47  double Original_Et = originalTower->et();
48  double etnew = Original_Et - (*emean_.find(-it)).second - (*esigma_.find(-it)).second;
49  if (etnew < 0.)
50  etnew = 0;
51  newjetet = newjetet + etnew;
52  jetOffset_[ijet] += Original_Et - etnew;
53  }
54  }
55 }
int def(FILE *, FILE *, int)
std::vector< double > jetOffset_
std::map< int, double > esigma_
std::vector< fastjet::PseudoJet > * fjJets_
std::vector< fastjet::PseudoJet > fjOriginalInputs_
int ieta(const reco::CandidatePtr &in) const
U second(std::pair< T, U > const &p)
std::vector< fastjet::PseudoJet > * fjInputs_
std::map< int, double > emean_
ActiveAreaSpecPtr fjActiveArea_
ClusterSequencePtr fjClusterSeq_
void rescaleRMS(double s)
void subtractPedestal(std::vector< fastjet::PseudoJet > &coll) override
#define LogDebug(id)
void ReflectedIterator::rescaleRMS ( double  s)

Definition at line 8 of file ReflectedIterator.cc.

8  {
9  for (std::map<int, double>::iterator iter = esigma_.begin(); iter != esigma_.end(); ++iter) {
10  iter->second = s * (iter->second);
11  }
12 }
std::map< int, double > esigma_
void ReflectedIterator::subtractPedestal ( std::vector< fastjet::PseudoJet > &  coll)
overridevirtual

Reimplemented from PileUpSubtractor.

Definition at line 57 of file ReflectedIterator.cc.

References LogDebug.

57  {
58  LogDebug("PileUpSubtractor") << "The subtractor subtracting pedestals...\n";
59 
60  int it = -100;
61 
62  vector<fastjet::PseudoJet> newcoll;
63 
64  for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin(), fjInputsEnd = coll.end();
65  input_object != fjInputsEnd;
66  ++input_object) {
67  reco::CandidatePtr const& itow = (*inputs_)[input_object->user_index()];
68 
69  it = ieta(itow);
70  iphi(itow);
71 
72  double Original_Et = itow->et();
73  if (sumRecHits_) {
74  Original_Et = getEt(itow);
75  }
76 
77  double etnew = Original_Et - (*(emean_.find(-it))).second - (*(esigma_.find(-it))).second;
78  float mScale = etnew / input_object->Et();
79  if (etnew < 0.)
80  mScale = 0.;
81 
82  math::XYZTLorentzVectorD towP4(input_object->px() * mScale,
83  input_object->py() * mScale,
84  input_object->pz() * mScale,
85  input_object->e() * mScale);
86 
87  int index = input_object->user_index();
88  input_object->reset(towP4.px(), towP4.py(), towP4.pz(), towP4.energy());
89  input_object->set_user_index(index);
90 
91  if (etnew > 0. && dropZeroTowers_)
92  newcoll.push_back(*input_object);
93  }
94 
95  if (dropZeroTowers_)
96  coll = newcoll;
97 }
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_
int iphi(const reco::CandidatePtr &in) const
double getEt(const reco::CandidatePtr &in) const
#define LogDebug(id)

Member Data Documentation

bool ReflectedIterator::dropZeroTowers_

Definition at line 22 of file ReflectedIterator.h.

bool ReflectedIterator::sumRecHits_

Definition at line 21 of file ReflectedIterator.h.