CMS 3D CMS Logo

MultipleAlgoIterator.cc
Go to the documentation of this file.
1 #include <memory>
2 
7 
9 using namespace std;
10 
12  : PileUpSubtractor(iConfig, std::move(iC)),
13  minimumTowersFraction_(iConfig.getParameter<double>("minimumTowersFraction")),
14  sumRecHits_(iConfig.getParameter<bool>("sumRecHits")),
15  dropZeroTowers_(iConfig.getUntrackedParameter<bool>("dropZeroTowers", true)) {
16  LogDebug("PileUpSubtractor") << "LIMITING THE MINIMUM TOWERS FRACTION TO : " << minimumTowersFraction_ << endl;
17 }
18 
20  for (std::map<int, double>::iterator iter = esigma_.begin(); iter != esigma_.end(); ++iter) {
21  iter->second = s * (iter->second);
22  }
23 }
24 
26  LogDebug("PileUpSubtractor") << "The subtractor correcting jets...\n";
27  jetOffset_.clear();
28 
29  using namespace reco;
30 
31  (*fjInputs_) = fjOriginalInputs_;
34  const fastjet::JetDefinition& def = *fjJetDefinition_;
35  if (!doAreaFastjet_ && !doRhoFastjet_) {
36  fjClusterSeq_ = std::make_shared<fastjet::ClusterSequence>(*fjInputs_, def);
37  } else {
38  fjClusterSeq_ = ClusterSequencePtr(new fastjet::ClusterSequenceArea(*fjInputs_, def, *fjActiveArea_));
39  }
40 
41  (*fjJets_) = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
42 
43  jetOffset_.reserve(fjJets_->size());
44 
45  vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin(), jetsEnd = fjJets_->end();
46  for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
47  int ijet = pseudojetTMP - fjJets_->begin();
48  jetOffset_[ijet] = 0;
49 
50  std::vector<fastjet::PseudoJet> towers = sorted_by_pt(pseudojetTMP->constituents());
51 
52  double newjetet = 0.;
53  for (vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(), towEnd = towers.end(); ito != towEnd; ++ito) {
54  const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
55  int it = ieta(originalTower);
56  double Original_Et = originalTower->et();
57  double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second;
58  if (etnew < 0.)
59  etnew = 0;
60  newjetet = newjetet + etnew;
61  jetOffset_[ijet] += Original_Et - etnew;
62  }
63  }
64 }
65 
66 void MultipleAlgoIterator::subtractPedestal(vector<fastjet::PseudoJet>& coll) {
67  LogDebug("PileUpSubtractor") << "The subtractor subtracting pedestals...\n";
68 
69  int it = -100;
70 
71  vector<fastjet::PseudoJet> newcoll;
72 
73  for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin(), fjInputsEnd = coll.end();
74  input_object != fjInputsEnd;
75  ++input_object) {
76  reco::CandidatePtr const& itow = (*inputs_)[input_object->user_index()];
77 
78  it = ieta(itow);
79  iphi(itow);
80 
81  double Original_Et = itow->et();
82  if (sumRecHits_) {
83  Original_Et = getEt(itow);
84  }
85 
86  double etnew = Original_Et - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
87  float mScale = etnew / input_object->Et();
88  if (etnew < 0.)
89  mScale = 0.;
90 
91  math::XYZTLorentzVectorD towP4(input_object->px() * mScale,
92  input_object->py() * mScale,
93  input_object->pz() * mScale,
94  input_object->e() * mScale);
95 
96  int index = input_object->user_index();
97  input_object->reset_momentum(towP4.px(), towP4.py(), towP4.pz(), towP4.energy());
98  input_object->set_user_index(index);
99 
100  if (etnew > 0. && dropZeroTowers_)
101  newcoll.push_back(*input_object);
102  }
103 
104  if (dropZeroTowers_)
105  coll = newcoll;
106 }
107 
108 void MultipleAlgoIterator::calculatePedestal(vector<fastjet::PseudoJet> const& coll) {
109  LogDebug("PileUpSubtractor") << "The subtractor calculating pedestals...\n";
110 
111  map<int, double> emean2;
112  map<int, int> ntowers;
113 
114  int ietaold = -10000;
115  int ieta0 = -100;
116 
117  // Initial values for emean_, emean2, esigma_, ntowers
118 
119  for (int i = ietamin_; i < ietamax_ + 1; i++) {
120  emean_[i] = 0.;
121  emean2[i] = 0.;
122  esigma_[i] = 0.;
123  ntowers[i] = 0;
124  }
125 
126  for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin(), fjInputsEnd = coll.end();
127  input_object != fjInputsEnd;
128  ++input_object) {
129  const reco::CandidatePtr& originalTower = (*inputs_)[input_object->user_index()];
130  ieta0 = ieta(originalTower);
131  double Original_Et = originalTower->et();
132  if (sumRecHits_) {
133  Original_Et = getEt(originalTower);
134  }
135 
136  if (ieta0 - ietaold != 0) {
137  emean_[ieta0] = emean_[ieta0] + Original_Et;
138  emean2[ieta0] = emean2[ieta0] + Original_Et * Original_Et;
139  ntowers[ieta0] = 1;
140  ietaold = ieta0;
141  } else {
142  emean_[ieta0] = emean_[ieta0] + Original_Et;
143  emean2[ieta0] = emean2[ieta0] + Original_Et * Original_Et;
144  ntowers[ieta0]++;
145  }
146  }
147 
148  for (map<int, int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++) {
149  int it = (*gt).first;
150 
151  double e1 = (*(emean_.find(it))).second;
152  double e2 = (*emean2.find(it)).second;
153  int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
154 
155  LogDebug("PileUpSubtractor") << " ieta : " << it << " number of towers : " << nt << " e1 : " << e1 << " e2 : " << e2
156  << "\n";
157 
158  if (nt > 0) {
159  emean_[it] = e1 / nt;
160  double eee = e2 / nt - e1 * e1 / (nt * nt);
161  if (eee < 0.)
162  eee = 0.;
163  esigma_[it] = nSigmaPU_ * sqrt(eee);
164  } else {
165  emean_[it] = 0.;
166  esigma_[it] = 0.;
167  }
168  LogDebug("PileUpSubtractor") << " ieta : " << it << " Pedestals : " << emean_[it] << " " << esigma_[it] << "\n";
169  }
170 }
171 
172 void MultipleAlgoIterator::calculateOrphanInput(vector<fastjet::PseudoJet>& orphanInput) {
173  LogDebug("PileUpSubtractor") << "The subtractor calculating orphan input...\n";
174 
175  (*fjInputs_) = fjOriginalInputs_;
176 
177  vector<int> jettowers; // vector of towers indexed by "user_index"
178  vector<pair<int, int> > excludedTowers; // vector of excluded ieta, iphi values
179 
180  for (auto const& pseudojetTMP : *fjJets_) {
181  if (pseudojetTMP.perp() < puPtMin_)
182  continue;
183 
184  // find towers within radiusPU_ of this jet
185  for (auto const im : allgeomid_) {
186  double dr = reco::deltaR(geo_->getPosition(im), pseudojetTMP);
187  auto exclude = find(excludedTowers.begin(), excludedTowers.end(), pair<int, int>(im.ieta(), im.iphi()));
188  if (dr < radiusPU_ && exclude == excludedTowers.end() &&
189  (geomtowers_[im.ieta()] - ntowersWithJets_[im.ieta()]) > minimumTowersFraction_ * (geomtowers_[im.ieta()])) {
190  ntowersWithJets_[im.ieta()]++;
191 
192  excludedTowers.push_back(pair<int, int>(im.ieta(), im.iphi()));
193  }
194  }
195 
196  for (auto const& it : *fjInputs_) {
197  int index = it.user_index();
198  int ie = ieta((*inputs_)[index]);
199  int ip = iphi((*inputs_)[index]);
200  auto exclude = find(excludedTowers.begin(), excludedTowers.end(), pair<int, int>(ie, ip));
201  if (exclude != excludedTowers.end()) {
202  jettowers.push_back(index);
203  }
204  } // initial input collection
205 
206  } // pseudojets
207 
208  //
209  // Create a new collections from the towers not included in jets
210  //
211 
212  for (auto const& it : *fjInputs_) {
213  int index = it.user_index();
214  vector<int>::const_iterator itjet = find(jettowers.begin(), jettowers.end(), index);
215  if (itjet == jettowers.end()) {
216  const reco::CandidatePtr& originalTower = (*inputs_)[index];
217  fastjet::PseudoJet orphan(originalTower->px(), originalTower->py(), originalTower->pz(), originalTower->energy());
218  orphan.set_user_index(index);
219 
220  orphanInput.push_back(orphan);
221  }
222  }
223 }
224 
226  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
227  const GlobalPoint& pos = geo_->getPosition(ctc->id());
228  double energy = ctc->emEnergy() + ctc->hadEnergy();
229  double et = energy * sin(pos.theta());
230  return et;
231 }
232 
234  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
235  const GlobalPoint& pos = geo_->getPosition(ctc->id());
236  double eta = pos.eta();
237  return eta;
238 }
PileUpSubtractor::ClusterSequencePtr
std::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
Definition: PileUpSubtractor.h:25
reco::Candidate::energy
virtual double energy() const =0
energy
CaloTower::emEnergy
double emEnergy() const
Definition: CaloTower.h:130
HLT_FULL_cff.towers
towers
Definition: HLT_FULL_cff.py:36379
PileUpSubtractor::ieta
int ieta(const reco::CandidatePtr &in) const
Definition: PileUpSubtractor.cc:325
PileUpSubtractor::iphi
int iphi(const reco::CandidatePtr &in) const
Definition: PileUpSubtractor.cc:336
electrons_cff.bool
bool
Definition: electrons_cff.py:366
mps_fire.i
i
Definition: mps_fire.py:428
MultipleAlgoIterator::calculatePedestal
void calculatePedestal(std::vector< fastjet::PseudoJet > const &coll) override
Definition: MultipleAlgoIterator.cc:108
MessageLogger.h
nt
int nt
Definition: AMPTWrapper.h:42
CaloGeometry::getPosition
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
PileUpSubtractor::puPtMin_
double puPtMin_
Definition: PileUpSubtractor.h:67
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
CaloTower::id
CaloTowerDetId id() const
Definition: CaloTower.h:123
PileUpSubtractor::doRhoFastjet_
bool doRhoFastjet_
Definition: PileUpSubtractor.h:65
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:222
PileUpSubtractor::allgeomid_
std::vector< HcalDetId > allgeomid_
Definition: PileUpSubtractor.h:79
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:46
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
CandidateFwd.h
PileUpSubtractor::fjActiveArea_
ActiveAreaSpecPtr fjActiveArea_
Definition: PileUpSubtractor.h:75
MultipleAlgoIterator::offsetCorrectJets
void offsetCorrectJets() override
Definition: MultipleAlgoIterator.cc:25
reco::Candidate::pz
virtual double pz() const =0
z coordinate of momentum vector
alignCSCRings.s
s
Definition: alignCSCRings.py:92
submitPVValidationJobs.gt
list gt
Definition: submitPVValidationJobs.py:663
MultipleAlgoIterator::getEta
double getEta(const reco::CandidatePtr &in) const
Definition: MultipleAlgoIterator.cc:233
MultipleAlgoIterator::rescaleRMS
void rescaleRMS(double s)
Definition: MultipleAlgoIterator.cc:19
PVValHelper::eta
Definition: PVValidationHelpers.h:70
PileUpSubtractor::doAreaFastjet_
bool doAreaFastjet_
Definition: PileUpSubtractor.h:64
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
CaloTower::hadEnergy
double hadEnergy() const
Definition: CaloTower.h:131
MultipleAlgoIterator::MultipleAlgoIterator
MultipleAlgoIterator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
Definition: MultipleAlgoIterator.cc:11
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
Point3DBase< float, GlobalTag >
PileUpSubtractor::geo_
CaloGeometry const * geo_
Definition: PileUpSubtractor.h:76
MultipleAlgoIterator::calculateOrphanInput
void calculateOrphanInput(std::vector< fastjet::PseudoJet > &orphanInput) override
Definition: MultipleAlgoIterator.cc:172
reco::Candidate::py
virtual double py() const =0
y coordinate of momentum vector
PileUpSubtractor::jetOffset_
std::vector< double > jetOffset_
Definition: PileUpSubtractor.h:85
funct::true
true
Definition: Factorize.h:173
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
PileUpSubtractor::fjOriginalInputs_
std::vector< fastjet::PseudoJet > fjOriginalInputs_
Definition: PileUpSubtractor.h:60
edm::ParameterSet
Definition: ParameterSet.h:47
MultipleAlgoIterator.h
deltaR.h
PileUpSubtractor
Definition: PileUpSubtractor.h:23
PileUpSubtractor::fjJets_
std::vector< fastjet::PseudoJet > * fjJets_
Definition: PileUpSubtractor.h:59
recoMuon::in
Definition: RecoMuonEnumerators.h:6
PileUpSubtractor::fjJetDefinition_
JetDefPtr fjJetDefinition_
Definition: PileUpSubtractor.h:55
StorageManager_cfg.e1
e1
Definition: StorageManager_cfg.py:16
PileUpSubtractor::geomtowers_
std::map< int, int > geomtowers_
Definition: PileUpSubtractor.h:80
CaloTower
Definition: CaloTower.h:26
PileUpSubtractor::esigma_
std::map< int, double > esigma_
Definition: PileUpSubtractor.h:82
EgHLTOffHistBins_cfi.et
et
Definition: EgHLTOffHistBins_cfi.py:8
MultipleAlgoIterator::getEt
double getEt(const reco::CandidatePtr &in) const
Definition: MultipleAlgoIterator.cc:225
MultipleAlgoIterator::dropZeroTowers_
bool dropZeroTowers_
Definition: MultipleAlgoIterator.h:20
MultipleAlgoIterator::sumRecHits_
bool sumRecHits_
Definition: MultipleAlgoIterator.h:19
edm::Ptr< Candidate >
PileUpSubtractor::fjInputs_
std::vector< fastjet::PseudoJet > * fjInputs_
Definition: PileUpSubtractor.h:58
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
PileUpSubtractor::ietamin_
int ietamin_
Definition: PileUpSubtractor.h:78
PileUpSubtractor::ietamax_
int ietamax_
Definition: PileUpSubtractor.h:77
reco::Candidate::px
virtual double px() const =0
x coordinate of momentum vector
PileUpSubtractor::nSigmaPU_
double nSigmaPU_
Definition: PileUpSubtractor.h:73
MultipleAlgoIterator::subtractPedestal
void subtractPedestal(std::vector< fastjet::PseudoJet > &coll) override
Definition: MultipleAlgoIterator.cc:66
PileUpSubtractor::inputs_
std::vector< edm::Ptr< reco::Candidate > > * inputs_
Definition: PileUpSubtractor.h:57
flavorHistoryFilter_cfi.dr
dr
Definition: flavorHistoryFilter_cfi.py:37
PileUpSubtractor::jetPtMin_
double jetPtMin_
Definition: PileUpSubtractor.h:66
reco::deltaR
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
spu::def
int def(FILE *, FILE *, int)
Definition: SherpackUtilities.cc:14
Candidate.h
MultipleAlgoIterator::minimumTowersFraction_
double minimumTowersFraction_
Definition: MultipleAlgoIterator.h:18
PileUpSubtractor::ntowersWithJets_
std::map< int, int > ntowersWithJets_
Definition: PileUpSubtractor.h:81
PileUpSubtractor::radiusPU_
double radiusPU_
Definition: PileUpSubtractor.h:74
PileUpSubtractor::emean_
std::map< int, double > emean_
Definition: PileUpSubtractor.h:83
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
PileUpSubtractor::fjClusterSeq_
ClusterSequencePtr fjClusterSeq_
Definition: PileUpSubtractor.h:56
reco::Candidate::et
virtual double et() const =0
transverse energy