CMS 3D CMS Logo

ReflectedIterator.cc
Go to the documentation of this file.
4 
6 using namespace std;
7 
9  for (std::map<int, double>::iterator iter = esigma_.begin(); iter != esigma_.end(); ++iter) {
10  iter->second = s * (iter->second);
11  }
12 }
13 
15  LogDebug("PileUpSubtractor") << "The subtractor correcting jets...\n";
16  jetOffset_.clear();
17 
18  using namespace reco;
19 
20  (*fjInputs_) = fjOriginalInputs_;
21  rescaleRMS(nSigmaPU_);
22  subtractPedestal(*fjInputs_);
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 }
56 
57 void ReflectedIterator::subtractPedestal(vector<fastjet::PseudoJet>& coll) {
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 }
98 
99 void ReflectedIterator::calculatePedestal(vector<fastjet::PseudoJet> const& coll) {
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 }
162 
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 }
170 
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 }
CaloTower::emEnergy
double emEnergy() const
Definition: CaloTower.h:130
HLT_FULL_cff.towers
towers
Definition: HLT_FULL_cff.py:36428
ReflectedIterator::offsetCorrectJets
void offsetCorrectJets() override
Definition: ReflectedIterator.cc:14
mps_fire.i
i
Definition: mps_fire.py:428
MessageLogger.h
ReflectedIterator.h
nt
int nt
Definition: AMPTWrapper.h:42
ReflectedIterator::getEta
double getEta(const reco::CandidatePtr &in) const
Definition: ReflectedIterator.cc:171
math::XYZTLorentzVectorD
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
pos
Definition: PixelAliasList.h:18
ReflectedIterator::calculatePedestal
void calculatePedestal(std::vector< fastjet::PseudoJet > const &coll) override
Definition: ReflectedIterator.cc:99
CaloTower::id
CaloTowerDetId id() const
Definition: CaloTower.h:123
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:222
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
LEDCalibrationChannels.iphi
iphi
Definition: LEDCalibrationChannels.py:64
CandidateFwd.h
alignCSCRings.s
s
Definition: alignCSCRings.py:92
submitPVValidationJobs.gt
list gt
Definition: submitPVValidationJobs.py:663
PVValHelper::eta
Definition: PVValidationHelpers.h:69
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
CaloTower::hadEnergy
double hadEnergy() const
Definition: CaloTower.h:131
ReflectedIterator::getEt
double getEt(const reco::CandidatePtr &in) const
Definition: ReflectedIterator.cc:163
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
Point3DBase< float, GlobalTag >
LEDCalibrationChannels.ieta
ieta
Definition: LEDCalibrationChannels.py:63
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:223
recoMuon::in
Definition: RecoMuonEnumerators.h:6
StorageManager_cfg.e1
e1
Definition: StorageManager_cfg.py:16
CaloTower
Definition: CaloTower.h:26
EgHLTOffHistBins_cfi.et
et
Definition: EgHLTOffHistBins_cfi.py:8
edm::Ptr< Candidate >
std
Definition: JetResolutionObject.h:76
ReflectedIterator::rescaleRMS
void rescaleRMS(double s)
Definition: ReflectedIterator.cc:8
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
ReflectedIterator::subtractPedestal
void subtractPedestal(std::vector< fastjet::PseudoJet > &coll) override
Definition: ReflectedIterator.cc:57
spu::def
int def(FILE *, FILE *, int)
Definition: SherpackUtilities.cc:14
Candidate.h
reco::Candidate::et
virtual double et() const =0
transverse energy