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