CMS 3D CMS Logo

MultipleAlgoIterator.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 = *fjJetDefinition_;
24  if (!doAreaFastjet_ && !doRhoFastjet_) {
25  fjClusterSeq_ = ClusterSequencePtr(new fastjet::ClusterSequence(*fjInputs_, def));
26  } else {
27  fjClusterSeq_ = ClusterSequencePtr(new fastjet::ClusterSequenceArea(*fjInputs_, def, *fjActiveArea_));
28  }
29 
30  (*fjJets_) = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
31 
32  jetOffset_.reserve(fjJets_->size());
33 
34  vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin(), jetsEnd = fjJets_->end();
35  for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
36  int ijet = pseudojetTMP - fjJets_->begin();
37  jetOffset_[ijet] = 0;
38 
39  std::vector<fastjet::PseudoJet> towers = sorted_by_pt(pseudojetTMP->constituents());
40 
41  double newjetet = 0.;
42  for (vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(), towEnd = towers.end(); ito != towEnd; ++ito) {
43  const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
44  int it = ieta(originalTower);
45  double Original_Et = originalTower->et();
46  double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second;
47  if (etnew < 0.)
48  etnew = 0;
49  newjetet = newjetet + etnew;
50  jetOffset_[ijet] += Original_Et - etnew;
51  }
52  }
53 }
54 
55 void MultipleAlgoIterator::subtractPedestal(vector<fastjet::PseudoJet>& coll) {
56  LogDebug("PileUpSubtractor") << "The subtractor subtracting pedestals...\n";
57 
58  int it = -100;
59 
60  vector<fastjet::PseudoJet> newcoll;
61 
62  for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin(), fjInputsEnd = coll.end();
63  input_object != fjInputsEnd;
64  ++input_object) {
65  reco::CandidatePtr const& itow = (*inputs_)[input_object->user_index()];
66 
67  it = ieta(itow);
68  iphi(itow);
69 
70  double Original_Et = itow->et();
71  if (sumRecHits_) {
72  Original_Et = getEt(itow);
73  }
74 
75  double etnew = Original_Et - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
76  float mScale = etnew / input_object->Et();
77  if (etnew < 0.)
78  mScale = 0.;
79 
80  math::XYZTLorentzVectorD towP4(input_object->px() * mScale,
81  input_object->py() * mScale,
82  input_object->pz() * mScale,
83  input_object->e() * mScale);
84 
85  int index = input_object->user_index();
86  input_object->reset_momentum(towP4.px(), towP4.py(), towP4.pz(), towP4.energy());
87  input_object->set_user_index(index);
88 
89  if (etnew > 0. && dropZeroTowers_)
90  newcoll.push_back(*input_object);
91  }
92 
93  if (dropZeroTowers_)
94  coll = newcoll;
95 }
96 
97 void MultipleAlgoIterator::calculatePedestal(vector<fastjet::PseudoJet> const& coll) {
98  LogDebug("PileUpSubtractor") << "The subtractor calculating pedestals...\n";
99 
100  map<int, double> emean2;
101  map<int, int> ntowers;
102 
103  int ietaold = -10000;
104  int ieta0 = -100;
105 
106  // Initial values for emean_, emean2, esigma_, ntowers
107 
108  for (int i = ietamin_; i < ietamax_ + 1; i++) {
109  emean_[i] = 0.;
110  emean2[i] = 0.;
111  esigma_[i] = 0.;
112  ntowers[i] = 0;
113  }
114 
115  for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin(), fjInputsEnd = coll.end();
116  input_object != fjInputsEnd;
117  ++input_object) {
118  const reco::CandidatePtr& originalTower = (*inputs_)[input_object->user_index()];
119  ieta0 = ieta(originalTower);
120  double Original_Et = originalTower->et();
121  if (sumRecHits_) {
122  Original_Et = getEt(originalTower);
123  }
124 
125  if (ieta0 - ietaold != 0) {
126  emean_[ieta0] = emean_[ieta0] + Original_Et;
127  emean2[ieta0] = emean2[ieta0] + Original_Et * Original_Et;
128  ntowers[ieta0] = 1;
129  ietaold = ieta0;
130  } else {
131  emean_[ieta0] = emean_[ieta0] + Original_Et;
132  emean2[ieta0] = emean2[ieta0] + Original_Et * Original_Et;
133  ntowers[ieta0]++;
134  }
135  }
136 
137  for (map<int, int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++) {
138  int it = (*gt).first;
139 
140  double e1 = (*(emean_.find(it))).second;
141  double e2 = (*emean2.find(it)).second;
142  int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
143 
144  LogDebug("PileUpSubtractor") << " ieta : " << it << " number of towers : " << nt << " e1 : " << e1 << " e2 : " << e2
145  << "\n";
146 
147  if (nt > 0) {
148  emean_[it] = e1 / nt;
149  double eee = e2 / nt - e1 * e1 / (nt * nt);
150  if (eee < 0.)
151  eee = 0.;
152  esigma_[it] = nSigmaPU_ * sqrt(eee);
153  } else {
154  emean_[it] = 0.;
155  esigma_[it] = 0.;
156  }
157  LogDebug("PileUpSubtractor") << " ieta : " << it << " Pedestals : " << emean_[it] << " " << esigma_[it] << "\n";
158  }
159 }
160 
162  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
163  const GlobalPoint& pos = geo_->getPosition(ctc->id());
164  double energy = ctc->emEnergy() + ctc->hadEnergy();
165  double et = energy * sin(pos.theta());
166  return et;
167 }
168 
170  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
171  const GlobalPoint& pos = geo_->getPosition(ctc->id());
172  double eta = pos.eta();
173  return eta;
174 }
#define LogDebug(id)
void calculatePedestal(std::vector< fastjet::PseudoJet > const &coll) override
double getEt(const reco::CandidatePtr &in) const
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
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
double getEta(const reco::CandidatePtr &in) const
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
U second(std::pair< T, U > const &p)
virtual double et() const =0
transverse energy
double emEnergy() const
Definition: CaloTower.h:134
T sqrt(T t)
Definition: SSEVec.h:19
void offsetCorrectJets() override
int nt
Definition: AMPTWrapper.h:42
void subtractPedestal(std::vector< fastjet::PseudoJet > &coll) override
double hadEnergy() const
Definition: CaloTower.h:135
CaloTowerDetId id() const
Definition: CaloTower.h:127
std::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
JetCorrectorParametersCollection coll
Definition: classes.h:10
T eta() const
Definition: PV3DBase.h:73
fixed size matrix
JetCorrectorParameters::Definitions def
Definition: classes.h:6