CMS 3D CMS Logo

PileUpSubtractor.cc
Go to the documentation of this file.
1 
3 
8 
14 
15 #include <map>
16 using namespace std;
17 
19 
20  geo_ = 0;
21  doAreaFastjet_ = iConfig.getParameter<bool>("doAreaFastjet");
22  doRhoFastjet_ = iConfig.getParameter<bool>("doRhoFastjet");
23  nSigmaPU_ = iConfig.getParameter<double>("nSigmaPU");
24  radiusPU_ = iConfig.getParameter<double>("radiusPU");
25  puPtMin_ = iConfig.getParameter<double>("puPtMin");
26  ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
27  activeAreaRepeats = iConfig.getParameter<int>("Active_Area_Repeats");
28  ghostArea = iConfig.getParameter<double>("GhostArea");
29 
30  if ( doAreaFastjet_ || doRhoFastjet_ ) {
31  fjActiveArea_ = ActiveAreaSpecPtr(new fastjet::ActiveAreaSpec(ghostEtaMax,
32  activeAreaRepeats,
33  ghostArea));
34  fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(ghostEtaMax) );
35  if ( ( ghostEtaMax < 0 ) || ( activeAreaRepeats < 0 ) || ( ghostArea < 0 ) )
36  throw cms::Exception("doAreaFastjet or doRhoFastjet") << "Parameters ghostEtaMax, activeAreaRepeats or ghostArea for doAreaFastjet/doRhoFastjet are not defined." << std::endl;
37  }
38 
39 }
40 
42  std::vector<fastjet::PseudoJet>& towers,
43  std::vector<fastjet::PseudoJet>& output){
44 
45  inputs_ = &input;
46  fjInputs_ = &towers;
47  fjJets_ = &output;
48  fjOriginalInputs_ = (*fjInputs_);
49  for (unsigned int i = 0; i < fjInputs_->size(); ++i){
50  fjOriginalInputs_[i].set_user_index((*fjInputs_)[i].user_index());
51  }
52 
53 }
54 
56  fjJetDefinition_ = JetDefPtr( new fastjet::JetDefinition( *jetDef ) );
57 }
58 
60 {
61 
62  LogDebug("PileUpSubtractor")<<"The subtractor setting up geometry...\n";
63 
64  if(geo_ == 0) {
66  iSetup.get<CaloGeometryRecord>().get(pG);
67  geo_ = pG.product();
68  std::vector<DetId> alldid = geo_->getValidDetIds();
69 
70  int ietaold = -10000;
71  ietamax_ = -10000;
72  ietamin_ = 10000;
73  for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
74  if( (*did).det() == DetId::Hcal ){
75  HcalDetId hid = HcalDetId(*did);
76  if( (hid).depth() == 1 ) {
77  allgeomid_.push_back(*did);
78 
79  if((hid).ieta() != ietaold){
80  ietaold = (hid).ieta();
81  geomtowers_[(hid).ieta()] = 1;
82  if((hid).ieta() > ietamax_) ietamax_ = (hid).ieta();
83  if((hid).ieta() < ietamin_) ietamin_ = (hid).ieta();
84  }
85  else{
86  geomtowers_[(hid).ieta()]++;
87  }
88  }
89  }
90  }
91  }
92 
93  for (int i = ietamin_; i<ietamax_+1; i++) {
94  ntowersWithJets_[i] = 0;
95  }
96 }
97 
98 //
99 // Calculate mean E and sigma from jet collection "coll".
100 //
101 void PileUpSubtractor::calculatePedestal( vector<fastjet::PseudoJet> const & coll )
102 {
103  LogDebug("PileUpSubtractor")<<"The subtractor calculating pedestals...\n";
104  map<int,double> emean2;
105  map<int,int> ntowers;
106 
107  int ietaold = -10000;
108  int ieta0 = -100;
109 
110  // Initial values for emean_, emean2, esigma_, ntowers
111 
112  for(int i = ietamin_; i < ietamax_+1; i++)
113  {
114  emean_[i] = 0.;
115  emean2[i] = 0.;
116  esigma_[i] = 0.;
117  ntowers[i] = 0;
118  }
119 
120  for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin (),
121  fjInputsEnd = coll.end();
122  input_object != fjInputsEnd; ++input_object) {
123  const reco::CandidatePtr & originalTower=(*inputs_)[ input_object->user_index()];
124  ieta0 = ieta( originalTower );
125  double Original_Et = originalTower->et();
126  if( ieta0-ietaold != 0 )
127  {
128  emean_[ieta0] = emean_[ieta0]+Original_Et;
129  emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
130  ntowers[ieta0] = 1;
131  ietaold = ieta0;
132  }
133  else
134  {
135  emean_[ieta0] = emean_[ieta0]+Original_Et;
136  emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
137  ntowers[ieta0]++;
138  }
139  }
140 
141  for(map<int,int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++)
142  {
143  int it = (*gt).first;
144 
145  double e1 = (*(emean_.find(it))).second;
146  double e2 = (*emean2.find(it)).second;
147  int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
148 
149  LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" number of towers : "<<nt<<" e1 : "<<e1<<" e2 : "<<e2<<"\n";
150  if(nt > 0) {
151  emean_[it] = e1/nt;
152  double eee = e2/nt - e1*e1/(nt*nt);
153  if(eee<0.) eee = 0.;
154  esigma_[it] = nSigmaPU_*sqrt(eee);
155  }
156  else
157  {
158  emean_[it] = 0.;
159  esigma_[it] = 0.;
160  }
161  LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" Pedestals : "<<emean_[it]<<" "<<esigma_[it]<<"\n";
162  }
163 }
164 
165 
166 //
167 // Subtract mean and sigma from fjInputs_
168 //
169 void PileUpSubtractor::subtractPedestal(vector<fastjet::PseudoJet> & coll)
170 {
171 
172  LogDebug("PileUpSubtractor")<<"The subtractor subtracting pedestals...\n";
173 
174  int it = -100;
175  for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
176  fjInputsEnd = coll.end();
177  input_object != fjInputsEnd; ++input_object) {
178 
179  reco::CandidatePtr const & itow = (*inputs_)[ input_object->user_index() ];
180 
181  it = ieta( itow );
182 
183  double etnew = itow->et() - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
184  float mScale = etnew/input_object->Et();
185  if(etnew < 0.) mScale = 0.;
186 
187  math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
188  input_object->pz()*mScale, input_object->e()*mScale);
189 
190  int index = input_object->user_index();
191  input_object->reset_momentum ( towP4.px(),
192  towP4.py(),
193  towP4.pz(),
194  towP4.energy() );
195  input_object->set_user_index(index);
196  }
197 }
198 
199 void PileUpSubtractor::calculateOrphanInput(vector<fastjet::PseudoJet> & orphanInput)
200 {
201 
202  LogDebug("PileUpSubtractor")<<"The subtractor calculating orphan input...\n";
203 
204  (*fjInputs_) = fjOriginalInputs_;
205 
206  vector<int> jettowers; // vector of towers indexed by "user_index"
207  vector<pair<int,int> > excludedTowers; // vector of excluded ieta, iphi values
208 
209  vector <fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
210  fjJetsEnd = fjJets_->end();
211  for (; pseudojetTMP != fjJetsEnd ; ++pseudojetTMP) {
212  if(pseudojetTMP->perp() < puPtMin_) continue;
213 
214  // find towers within radiusPU_ of this jet
215  for(vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++)
216  {
217  double dr = reco::deltaR(geo_->getPosition((DetId)(*im)),(*pseudojetTMP));
218  vector<pair<int,int> >::const_iterator exclude = find(excludedTowers.begin(),excludedTowers.end(),pair<int,int>(im->ieta(),im->iphi()));
219  if( dr < radiusPU_ && exclude == excludedTowers.end()) {
220  ntowersWithJets_[(*im).ieta()]++;
221  excludedTowers.push_back(pair<int,int>(im->ieta(),im->iphi()));
222  }
223  }
224  vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
225  fjInputsEnd = fjInputs_->end();
226 
227  for (; it != fjInputsEnd; ++it ) {
228  int index = it->user_index();
229  int ie = ieta((*inputs_)[index]);
230  int ip = iphi((*inputs_)[index]);
231  vector<pair<int,int> >::const_iterator exclude = find(excludedTowers.begin(),excludedTowers.end(),pair<int,int>(ie,ip));
232  if(exclude != excludedTowers.end()) {
233  jettowers.push_back(index);
234  } //dr < radiusPU_
235  } // initial input collection
236  } // pseudojets
237 
238  //
239  // Create a new collections from the towers not included in jets
240  //
241  for(vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
242  fjInputsEnd = fjInputs_->end(); it != fjInputsEnd; ++it ) {
243  int index = it->user_index();
244  vector<int>::const_iterator itjet = find(jettowers.begin(),jettowers.end(),index);
245  if( itjet == jettowers.end() ){
246  const reco::CandidatePtr& originalTower = (*inputs_)[index];
247  fastjet::PseudoJet orphan(originalTower->px(),originalTower->py(),originalTower->pz(),originalTower->energy());
248  orphan.set_user_index(index);
249 
250  orphanInput.push_back(orphan);
251  }
252  }
253 }
254 
255 
257 {
258  LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
259  jetOffset_.clear();
260  using namespace reco;
261 
262  //
263  // Reestimate energy of jet (energy of jet with initial map)
264  //
265  jetOffset_.reserve(fjJets_->size());
266  vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
267  jetsEnd = fjJets_->end();
268  for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
269  int ijet = pseudojetTMP - fjJets_->begin();
270  jetOffset_[ijet] = 0;
271 
272  std::vector<fastjet::PseudoJet> towers =
273  fastjet::sorted_by_pt( pseudojetTMP->constituents() );
274  double newjetet = 0.;
275  for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
276  towEnd = towers.end();
277  ito != towEnd;
278  ++ito)
279  {
280  const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
281  int it = ieta( originalTower );
282  double Original_Et = originalTower->et();
283  double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second;
284  if(etnew < 0.) etnew = 0;
285  newjetet = newjetet + etnew;
286  jetOffset_[ijet] += Original_Et - etnew;
287  }
288  double mScale = newjetet/pseudojetTMP->Et();
289  LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() : "<<pseudojetTMP->Et()<<"\n";
290  LogDebug("PileUpSubtractor")<<"newjetet : "<<newjetet<<"\n";
291  LogDebug("PileUpSubtractor")<<"jetOffset_[ijet] : "<<jetOffset_[ijet]<<"\n";
292  LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() - jetOffset_[ijet] : "<<pseudojetTMP->Et() - jetOffset_[ijet]<<"\n";
293  LogDebug("PileUpSubtractor")<<"Scale is : "<<mScale<<"\n";
294  int cshist = pseudojetTMP->cluster_hist_index();
295  pseudojetTMP->reset_momentum(pseudojetTMP->px()*mScale, pseudojetTMP->py()*mScale,
296  pseudojetTMP->pz()*mScale, pseudojetTMP->e()*mScale);
297  pseudojetTMP->set_cluster_hist_index(cshist);
298 
299  }
300 }
301 
302 double PileUpSubtractor::getCone(double cone, double eta, double phi, double& et, double& pu){
303  pu = 0;
304 
305  for(vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++){
306  if( im->depth() != 1 ) continue;
307  const GlobalPoint& point = geo_->getPosition((DetId)(*im));
308  double dr = reco::deltaR(point.eta(),point.phi(),eta,phi);
309  if( dr < cone){
310  pu += (*emean_.find(im->ieta())).second+(*esigma_.find(im->ieta())).second;
311  }
312  }
313 
314  return pu;
315 }
316 
318  int it = ieta(in);
319  return (*emean_.find(it)).second;
320 }
321 
323  int it = ieta(in);
324  return (*esigma_.find(it)).second;
325 }
326 
328  int it = ieta(in);
329  return (*emean_.find(it)).second + (*esigma_.find(it)).second;
330 }
331 
333  int it = ieta(in);
334 
335  int n = (*(geomtowers_.find(it))).second - (*(ntowersWithJets_.find(it))).second;
336  return n;
337 
338 }
339 
341  int it = ieta(in);
342  int n = (*(ntowersWithJets_.find(it))).second;
343  return n;
344 
345 }
346 
347 
349  int it = 0;
350  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
351  if(ctc){
352  it = ctc->id().ieta();
353  } else
354  {
355  throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
356  }
357  return it;
358 }
359 
361  int it = 0;
362  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
363  if(ctc){
364  it = ctc->id().iphi();
365  } else
366  {
367  throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
368  }
369  return it;
370 }
371 
372 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
374 
376  desc.add<bool> ("doAreaFastjet", false);
377  desc.add<bool> ("doRhoFastjet", false);
378  desc.add<double> ("Ghost_EtaMax", 5);
379  desc.add<double> ("GhostArea", 0.01);
380  desc.add<int> ("Active_Area_Repeats", 1);
381  desc.add<double> ("puPtMin", 10.);
382  desc.add<double> ("nSigmaPU", 1.);
383  desc.add<double> ("radiusPU", 0.5);
384  descriptions.addDefault(desc);
385 }
386 
389 
390 
391 
#define LogDebug(id)
T getParameter(std::string const &) const
virtual void calculateOrphanInput(std::vector< fastjet::PseudoJet > &orphanInput)
boost::shared_ptr< fastjet::JetDefinition > JetDefPtr
virtual double getMeanAtTower(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
int getNwithJets(const reco::CandidatePtr &in) const
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:159
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual double getPileUpAtTower(const reco::CandidatePtr &in) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
static std::string const input
Definition: EdmProvDump.cc:44
int ieta(const reco::CandidatePtr &in) const
U second(std::pair< T, U > const &p)
virtual void offsetCorrectJets()
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
T sqrt(T t)
Definition: SSEVec.h:18
int iphi() const
get the tower iphi
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
ParameterDescriptionBase * add(U const &iLabel, T const &value)
boost::shared_ptr< fastjet::GhostedAreaSpec > ActiveAreaSpecPtr
virtual void subtractPedestal(std::vector< fastjet::PseudoJet > &coll)
int nt
Definition: AMPTWrapper.h:32
#define EDM_REGISTER_PLUGINFACTORY(_factory_, _category_)
Definition: PluginFactory.h:90
Float e1
Definition: deltaR.h:20
Definition: DetId.h:18
int iphi(const reco::CandidatePtr &in) const
CaloTowerDetId id() const
Definition: CaloTower.h:103
JetCorrectorParametersCollection coll
Definition: classes.h:10
virtual double getSigmaAtTower(const reco::CandidatePtr &in) const
const T & get() const
Definition: EventSetup.h:55
boost::shared_ptr< fastjet::RangeDefinition > RangeDefPtr
virtual void setDefinition(JetDefPtr const &jetDef)
et
define resolution functions of each parameter
Float e2
Definition: deltaR.h:21
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:92
T eta() const
Definition: PV3DBase.h:76
fixed size matrix
PileUpSubtractor(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
virtual void calculatePedestal(std::vector< fastjet::PseudoJet > const &coll)
virtual void setupGeometryMap(edm::Event &iEvent, const edm::EventSetup &iSetup)
int ieta() const
get the tower ieta
virtual double getCone(double cone, double eta, double phi, double &et, double &pu)
virtual void reset(std::vector< edm::Ptr< reco::Candidate > > &input, std::vector< fastjet::PseudoJet > &towers, std::vector< fastjet::PseudoJet > &output)
T const * product() const
Definition: ESHandle.h:86
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
int getN(const reco::CandidatePtr &in) const
boost::shared_ptr< fastjet::JetDefinition > JetDefPtr