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_ = nullptr;
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  jetPtMin_ = iConfig.getParameter<double>("jetPtMin");
26  puPtMin_ = iConfig.getParameter<double>("puPtMin");
27  ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
28  activeAreaRepeats = iConfig.getParameter<int>("Active_Area_Repeats");
29  ghostArea = iConfig.getParameter<double>("GhostArea");
30 
31  if ( doAreaFastjet_ || doRhoFastjet_ ) {
32  fjActiveArea_ = ActiveAreaSpecPtr(new fastjet::ActiveAreaSpec(ghostEtaMax,
33  activeAreaRepeats,
34  ghostArea));
35  fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(ghostEtaMax) );
36  if ( ( ghostEtaMax < 0 ) || ( activeAreaRepeats < 0 ) || ( ghostArea < 0 ) )
37  throw cms::Exception("doAreaFastjet or doRhoFastjet") << "Parameters ghostEtaMax, activeAreaRepeats or ghostArea for doAreaFastjet/doRhoFastjet are not defined." << std::endl;
38  }
39 
40 }
41 
43  std::vector<fastjet::PseudoJet>& towers,
44  std::vector<fastjet::PseudoJet>& output){
45 
46  inputs_ = &input;
47  fjInputs_ = &towers;
48  fjJets_ = &output;
49  fjOriginalInputs_ = (*fjInputs_);
50  for (unsigned int i = 0; i < fjInputs_->size(); ++i){
51  fjOriginalInputs_[i].set_user_index((*fjInputs_)[i].user_index());
52  }
53 
54 }
55 
57  fjJetDefinition_ = JetDefPtr( new fastjet::JetDefinition( *jetDef ) );
58 }
59 
61 {
62 
63  LogDebug("PileUpSubtractor")<<"The subtractor setting up geometry...\n";
64 
65  if(geo_ == nullptr) {
67  iSetup.get<CaloGeometryRecord>().get(pG);
68  geo_ = pG.product();
69  std::vector<DetId> alldid = geo_->getValidDetIds();
70 
71  int ietaold = -10000;
72  ietamax_ = -10000;
73  ietamin_ = 10000;
74  for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
75  if( (*did).det() == DetId::Hcal ){
76  HcalDetId hid = HcalDetId(*did);
77  if( (hid).depth() == 1 ) {
78  allgeomid_.push_back(*did);
79 
80  if((hid).ieta() != ietaold){
81  ietaold = (hid).ieta();
82  geomtowers_[(hid).ieta()] = 1;
83  if((hid).ieta() > ietamax_) ietamax_ = (hid).ieta();
84  if((hid).ieta() < ietamin_) ietamin_ = (hid).ieta();
85  }
86  else{
87  geomtowers_[(hid).ieta()]++;
88  }
89  }
90  }
91  }
92  }
93 
94  for (int i = ietamin_; i<ietamax_+1; i++) {
95  ntowersWithJets_[i] = 0;
96  }
97 }
98 
99 //
100 // Calculate mean E and sigma from jet collection "coll".
101 //
102 void PileUpSubtractor::calculatePedestal( vector<fastjet::PseudoJet> const & coll )
103 {
104  LogDebug("PileUpSubtractor")<<"The subtractor calculating pedestals...\n";
105  map<int,double> emean2;
106  map<int,int> ntowers;
107 
108  int ietaold = -10000;
109  int ieta0 = -100;
110 
111  // Initial values for emean_, emean2, esigma_, ntowers
112 
113  for(int i = ietamin_; i < ietamax_+1; i++)
114  {
115  emean_[i] = 0.;
116  emean2[i] = 0.;
117  esigma_[i] = 0.;
118  ntowers[i] = 0;
119  }
120 
121  for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin (),
122  fjInputsEnd = coll.end();
123  input_object != fjInputsEnd; ++input_object) {
124  const reco::CandidatePtr & originalTower=(*inputs_)[ input_object->user_index()];
125  ieta0 = ieta( originalTower );
126  double Original_Et = originalTower->et();
127  if( ieta0-ietaold != 0 )
128  {
129  emean_[ieta0] = emean_[ieta0]+Original_Et;
130  emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
131  ntowers[ieta0] = 1;
132  ietaold = ieta0;
133  }
134  else
135  {
136  emean_[ieta0] = emean_[ieta0]+Original_Et;
137  emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
138  ntowers[ieta0]++;
139  }
140  }
141 
142  for(map<int,int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++)
143  {
144  int it = (*gt).first;
145 
146  double e1 = (*(emean_.find(it))).second;
147  double e2 = (*emean2.find(it)).second;
148  int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
149 
150  LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" number of towers : "<<nt<<" e1 : "<<e1<<" e2 : "<<e2<<"\n";
151  if(nt > 0) {
152  emean_[it] = e1/nt;
153  double eee = e2/nt - e1*e1/(nt*nt);
154  if(eee<0.) eee = 0.;
155  esigma_[it] = nSigmaPU_*sqrt(eee);
156  }
157  else
158  {
159  emean_[it] = 0.;
160  esigma_[it] = 0.;
161  }
162  LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" Pedestals : "<<emean_[it]<<" "<<esigma_[it]<<"\n";
163  }
164 }
165 
166 
167 //
168 // Subtract mean and sigma from fjInputs_
169 //
170 void PileUpSubtractor::subtractPedestal(vector<fastjet::PseudoJet> & coll)
171 {
172 
173  LogDebug("PileUpSubtractor")<<"The subtractor subtracting pedestals...\n";
174 
175  int it = -100;
176  for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
177  fjInputsEnd = coll.end();
178  input_object != fjInputsEnd; ++input_object) {
179 
180  reco::CandidatePtr const & itow = (*inputs_)[ input_object->user_index() ];
181 
182  it = ieta( itow );
183 
184  double etnew = itow->et() - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
185  float mScale = etnew/input_object->Et();
186  if(etnew < 0.) mScale = 0.;
187 
188  math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
189  input_object->pz()*mScale, input_object->e()*mScale);
190 
191  int index = input_object->user_index();
192  input_object->reset_momentum ( towP4.px(),
193  towP4.py(),
194  towP4.pz(),
195  towP4.energy() );
196  input_object->set_user_index(index);
197  }
198 }
199 
200 void PileUpSubtractor::calculateOrphanInput(vector<fastjet::PseudoJet> & orphanInput)
201 {
202 
203  LogDebug("PileUpSubtractor")<<"The subtractor calculating orphan input...\n";
204 
205  (*fjInputs_) = fjOriginalInputs_;
206 
207  vector<int> jettowers; // vector of towers indexed by "user_index"
208  vector<pair<int,int> > excludedTowers; // vector of excluded ieta, iphi values
209 
210  vector <fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
211  fjJetsEnd = fjJets_->end();
212  for (; pseudojetTMP != fjJetsEnd ; ++pseudojetTMP) {
213  if(pseudojetTMP->perp() < puPtMin_) continue;
214 
215  // find towers within radiusPU_ of this jet
216  for(vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++)
217  {
218  double dr = reco::deltaR(geo_->getPosition((DetId)(*im)),(*pseudojetTMP));
219  vector<pair<int,int> >::const_iterator exclude = find(excludedTowers.begin(),excludedTowers.end(),pair<int,int>(im->ieta(),im->iphi()));
220  if( dr < radiusPU_ && exclude == excludedTowers.end()) {
221  ntowersWithJets_[(*im).ieta()]++;
222  excludedTowers.push_back(pair<int,int>(im->ieta(),im->iphi()));
223  }
224  }
225  vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
226  fjInputsEnd = fjInputs_->end();
227 
228  for (; it != fjInputsEnd; ++it ) {
229  int index = it->user_index();
230  int ie = ieta((*inputs_)[index]);
231  int ip = iphi((*inputs_)[index]);
232  vector<pair<int,int> >::const_iterator exclude = find(excludedTowers.begin(),excludedTowers.end(),pair<int,int>(ie,ip));
233  if(exclude != excludedTowers.end()) {
234  jettowers.push_back(index);
235  } //dr < radiusPU_
236  } // initial input collection
237  } // pseudojets
238 
239  //
240  // Create a new collections from the towers not included in jets
241  //
242  for(vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
243  fjInputsEnd = fjInputs_->end(); it != fjInputsEnd; ++it ) {
244  int index = it->user_index();
245  vector<int>::const_iterator itjet = find(jettowers.begin(),jettowers.end(),index);
246  if( itjet == jettowers.end() ){
247  const reco::CandidatePtr& originalTower = (*inputs_)[index];
248  fastjet::PseudoJet orphan(originalTower->px(),originalTower->py(),originalTower->pz(),originalTower->energy());
249  orphan.set_user_index(index);
250 
251  orphanInput.push_back(orphan);
252  }
253  }
254 }
255 
256 
258 {
259  LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
260  jetOffset_.clear();
261  using namespace reco;
262 
263  //
264  // Reestimate energy of jet (energy of jet with initial map)
265  //
266  jetOffset_.reserve(fjJets_->size());
267  vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
268  jetsEnd = fjJets_->end();
269  for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
270  int ijet = pseudojetTMP - fjJets_->begin();
271  jetOffset_[ijet] = 0;
272 
273  std::vector<fastjet::PseudoJet> towers =
274  fastjet::sorted_by_pt( pseudojetTMP->constituents() );
275  double newjetet = 0.;
276  for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
277  towEnd = towers.end();
278  ito != towEnd;
279  ++ito)
280  {
281  const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
282  int it = ieta( originalTower );
283  double Original_Et = originalTower->et();
284  double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second;
285  if(etnew < 0.) etnew = 0;
286  newjetet = newjetet + etnew;
287  jetOffset_[ijet] += Original_Et - etnew;
288  }
289  double mScale = newjetet/pseudojetTMP->Et();
290  LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() : "<<pseudojetTMP->Et()<<'\n';
291  LogDebug("PileUpSubtractor")<<"newjetet : "<<newjetet<<'\n';
292  LogDebug("PileUpSubtractor")<<"jetOffset_[ijet] : "<<jetOffset_[ijet]<<'\n';
293  LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() - jetOffset_[ijet] : "<<pseudojetTMP->Et() - jetOffset_[ijet]<<'\n';
294  LogDebug("PileUpSubtractor")<<"Scale is : "<<mScale<<'\n';
295  int cshist = pseudojetTMP->cluster_hist_index();
296  pseudojetTMP->reset_momentum(pseudojetTMP->px()*mScale, pseudojetTMP->py()*mScale,
297  pseudojetTMP->pz()*mScale, pseudojetTMP->e()*mScale);
298  pseudojetTMP->set_cluster_hist_index(cshist);
299 
300  }
301 }
302 
303 double PileUpSubtractor::getCone(double cone, double eta, double phi, double& et, double& pu){
304  pu = 0;
305 
306  for(vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++){
307  if( im->depth() != 1 ) continue;
308  const GlobalPoint& point = geo_->getPosition((DetId)(*im));
309  double dr = reco::deltaR(point.eta(),point.phi(),eta,phi);
310  if( dr < cone){
311  pu += (*emean_.find(im->ieta())).second+(*esigma_.find(im->ieta())).second;
312  }
313  }
314 
315  return pu;
316 }
317 
319  int it = ieta(in);
320  return (*emean_.find(it)).second;
321 }
322 
324  int it = ieta(in);
325  return (*esigma_.find(it)).second;
326 }
327 
329  int it = ieta(in);
330  return (*emean_.find(it)).second + (*esigma_.find(it)).second;
331 }
332 
334  int it = ieta(in);
335 
336  int n = (*(geomtowers_.find(it))).second - (*(ntowersWithJets_.find(it))).second;
337  return n;
338 
339 }
340 
342  int it = ieta(in);
343  int n = (*(ntowersWithJets_.find(it))).second;
344  return n;
345 
346 }
347 
348 
350  int it = 0;
351  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
352  if(ctc){
353  it = ctc->id().ieta();
354  } else
355  {
356  throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
357  }
358  return it;
359 }
360 
362  int it = 0;
363  const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
364  if(ctc){
365  it = ctc->id().iphi();
366  } else
367  {
368  throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
369  }
370  return it;
371 }
372 
375 
376 
377 
#define LogDebug(id)
std::shared_ptr< fastjet::JetDefinition > JetDefPtr
T getParameter(std::string const &) const
virtual void calculateOrphanInput(std::vector< fastjet::PseudoJet > &orphanInput)
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
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
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
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:91
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:58
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:97
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