CMS 3D CMS Logo

ParticleTowerProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ParticleTowerProducer
4 // Class: ParticleTowerProducer
5 //
13 //
14 // Original Author: Yetkin Yilmaz,32 4-A08,+41227673039,
15 // Created: Thu Jan 20 19:53:58 CET 2011
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 
24 // user include files
27 
30 
32 
34 
36 
41 
46 
47 #include "TMath.h"
48 #include "TRandom.h"
49 
50 
51 
52 
53 //
54 // constants, enums and typedefs
55 //
56 
57 
58 //
59 // static data member definitions
60 //
61 
62 //
63 // constructors and destructor
64 //
66  geo_(nullptr)
67 {
68  //register your products
69  src_ = consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("src"));
70  useHF_ = iConfig.getParameter<bool>("useHF");
71 
72  produces<CaloTowerCollection>();
73 
74  //now do what ever other initialization is needed
75  random_ = new TRandom();
76  PI = TMath::Pi();
77 
78 
79 
80 }
81 
82 
84 {
85 
86  // do anything here that needs to be done at desctruction time
87  // (e.g. close files, deallocate resources etc.)
88  delete random_;
89 }
90 
91 
92 //
93 // member functions
94 //
95 
96 // ------------ method called to produce the data ------------
97 void
99 {
100  using namespace edm;
101 
102  if(!geo_){
104  iSetup.get<CaloGeometryRecord>().get(pG);
105  geo_ = pG.product();
106  }
107 
108 
109  resetTowers(iEvent, iSetup);
110 
111 
113  iEvent.getByToken(src_, inputsHandle);
114 
115  for(reco::PFCandidateCollection::const_iterator ci = inputsHandle->begin(); ci!=inputsHandle->end(); ++ci) {
116 
117  const reco::PFCandidate& particle = *ci;
118 
119  // put a cutoff if you want
120  //if(particle.et() < 0.3) continue;
121 
122  double eta = particle.eta();
123 
124  if(!useHF_ && fabs(eta) > 3. ) continue;
125 
126  int ieta = eta2ieta(eta);
127  if(ieta==-1) continue;
128  if(eta<0) ieta *= -1;
129 
130  int iphi = phi2iphi(particle.phi(),ieta);
131 
133  if(abs(ieta)<=16) sd = HcalBarrel;
134  else if(fabs(eta)< 3.) sd = HcalEndcap; // Use the endcap until eta =3
135  else sd = HcalForward;
136 
137  HcalDetId hid = HcalDetId(sd,ieta,iphi,1); // assume depth=1
138 
139  // check against the old method (boundaries slightly shifted in the endcaps
140  /*
141  HcalDetId hid_orig = getNearestTower(particle);
142 
143  if(hid != hid_orig)
144  {
145  if(abs((hid).ieta()-(hid_orig).ieta())>1 || abs((hid).iphi()-(hid_orig).iphi()) >1 ){
146 
147  int phi_diff = 1;
148  if(abs((hid).ieta())>39 || abs((hid_orig).ieta())>39) phi_diff = 4;
149  else if(abs((hid).ieta())>20 || abs((hid_orig).ieta())>20) phi_diff = 2;
150 
151  if(abs((hid).ieta()-(hid_orig).ieta()) > 1 || abs((hid).iphi()-(hid_orig).iphi()) > phi_diff ){
152  std::cout<<" eta "<<eta<<std::endl;
153  std::cout<<" hid "<<hid<<std::endl;
154  std::cout<<" old method "<<hid_orig<<std::endl;
155  }
156  }
157  }
158  */
159  towers_[hid] += particle.et();
160  }
161 
162 
163  auto prod = std::make_unique<CaloTowerCollection>();
164 
165  for ( std::map< DetId, double >::const_iterator iter = towers_.begin();
166  iter != towers_.end(); ++iter ){
167 
168  HcalDetId id = HcalDetId(iter->first);
169  CaloTowerDetId newTowerId(id.ieta(),id.iphi());
170  double et = iter->second;
171 
172  if(et>0){
173 
174  GlobalPoint pos =geo_->getGeometry(newTowerId)->getPosition();
175 
176  if(!useHF_ && fabs(pos.eta()) > 3. ) continue;
177 
178  // currently sets et = pt, mass to zero
179  // pt, eta , phi, mass
180  reco::Particle::PolarLorentzVector p4(et,pos.eta(),pos.phi(),0.);
181 
182  CaloTower newTower(newTowerId,et,0,0,0,0,p4,pos,pos);
183  prod->push_back(newTower);
184  }
185  }
186 
187 
188  //For reference, Calo Tower Constructors
189 
190  /*
191  CaloTower(const CaloTowerDetId& id,
192  double emE, double hadE, double outerE,
193  int ecal_tp, int hcal_tp,
194  const PolarLorentzVector p4,
195  GlobalPoint emPosition, GlobalPoint hadPosition);
196 
197  CaloTower(const CaloTowerDetId& id,
198  double emE, double hadE, double outerE,
199  int ecal_tp, int hcal_tp,
200  const LorentzVector p4,
201  GlobalPoint emPosition, GlobalPoint hadPosition);
202  */
203 
204 
205  iEvent.put(std::move(prod));
206 
207 
208 }
209 
210 // ------------ method called once each job just before starting event loop ------------
211 void
213 {
214  // tower edges from fast sim, used starting at index 30 for the HF
215  const double etatow[42] = {0.000, 0.087, 0.174, 0.261, 0.348, 0.435, 0.522, 0.609, 0.696, 0.783, 0.870, 0.957, 1.044, 1.131, 1.218, 1.305, 1.392, 1.479, 1.566, 1.653, 1.740, 1.830, 1.930, 2.043, 2.172, 2.322, 2.500, 2.650, 2.853, 3.000, 3.139, 3.314, 3.489, 3.664, 3.839, 4.013, 4.191, 4.363, 4.538, 4.716, 4.889, 5.191};
216 
217  // tower centers derived directly from by looping over HCAL towers and printing out the values, used up until eta = 3 where the HF/Endcap interface becomes problematic
218  const double etacent[40] = {
219  0.0435, // 1
220  0.1305, // 2
221  0.2175, // 3
222  0.3045, // 4
223  0.3915, // 5
224  0.4785, // 6
225  0.5655, // 7
226  0.6525, // 8
227  0.7395, // 9
228  0.8265, // 10
229  0.9135, // 11
230  1.0005, // 12
231  1.0875, // 13
232  1.1745, // 14
233  1.2615, // 15
234  1.3485, // 16
235  1.4355, // 17
236  1.5225, // 18
237  1.6095, // 19
238  1.6965, // 20
239  1.785 , // 21
240  1.88 , // 22
241  1.9865, // 23
242  2.1075, // 24
243  2.247 , // 25
244  2.411 , // 26
245  2.575 , // 27
246  2.759 , // 28
247  2.934 , // 29
248  2.90138,// 29
249  3.04448,// 30
250  3.21932,// 31
251  3.39462,// 32
252  3.56966,// 33
253  3.74485,// 34
254  3.91956,// 35
255  4.09497,// 36
256  4.27007,// 37
257  4.44417,// 38
258  4.62046 // 39
259  };
260 
261  // Use the real towers centrers for the barrel and endcap up to eta=3
262  etaedge[0] = 0.;
263  for(int i=1;i<30;i++){
264  etaedge[i] = (etacent[i-1]-etaedge[i-1])*2.0 + etaedge[i-1];
265  //std::cout<<" i "<<i<<" etaedge "<<etaedge[i]<<std::endl;
266  }
267  // beyond eta = 3 just use the fast sim values
268  if(useHF_){
269  for(int i=30;i<42;i++){
270  etaedge[i]=etatow[i];
271  //std::cout<<" i "<<i<<" etaedge "<<etaedge[i]<<std::endl;
272  }
273  }
274 }
275 
276 // ------------ method called once each job just after ending the event loop ------------
277 void
279 }
280 
281 
283 {
284 
285  std::vector<DetId> alldid = geo_->getValidDetIds();
286 
287  for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
288  if( (*did).det() == DetId::Hcal ){
289  HcalDetId hid = HcalDetId(*did);
290  if( hid.depth() == 1 ) {
291 
292  if(!useHF_){
293  GlobalPoint pos =geo_->getGeometry(hid)->getPosition();
294  //if((hid).iphi()==1)std::cout<<" ieta "<<(hid).ieta()<<" eta "<<pos.eta()<<" iphi "<<(hid).iphi()<<" phi "<<pos.phi()<<std::endl;
295  if(fabs(pos.eta())>3.) continue;
296  }
297  towers_[(*did)] = 0.;
298  }
299 
300  }
301  }
302 
303 
304 
305 
306  //print out tower eta boundaries
307  /*
308  int current_ieta=0;
309  float testphi =0.;
310  for(double testeta=1.7; testeta<=5.4;testeta+=0.00001){
311  HcalDetId hid = getNearestTower(testeta,testphi);
312  if((hid).ieta()!=current_ieta) {
313  std::cout<<" new ieta "<<current_ieta<<" testeta "<<testeta<<std::endl;
314  current_ieta = (hid).ieta();
315  }
316  }
317  */
318 
319 }
320 
321 
322 
323 
325 
326  double eta = in.eta();
327  double phi = in.phi();
328 
329  double minDeltaR = 9999;
330 
331  std::vector<DetId> alldid = geo_->getValidDetIds();
332 
333  DetId returnId;
334 
335  //int nclosetowers=0;
336 
337  for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
338  if( (*did).det() == DetId::Hcal ){
339 
340  HcalDetId hid(*did);
341 
342  // which layer is irrelevant for an eta-phi map, no?
343 
344  if( hid.depth() != 1 ) continue;
345 
346  GlobalPoint pos =geo_->getGeometry(hid)->getPosition();
347 
348  double hcalEta = pos.eta();
349  double hcalPhi = pos.phi();
350 
351  //std::cout<<" transverse distance = "<<sqrt(pos.x()*pos.x()+pos.y()*pos.y())<<std::endl;
352 
353  //std::cout<<" ieta "<<(hid).ieta()<<" iphi "<<(hid).iphi()<<" hcalEta "<<hcalEta<<" hcalPhi "<<hcalPhi<<std::endl;
354 
355  double deltaR = reco::deltaR(eta,phi,hcalEta,hcalPhi);
356 
357  // need to factor in the size of the tower
358  double towersize = 0.087;
359 
360  int ieta = (hid).ieta();
361 
362  if(abs(ieta)>21){
363  if(abs(ieta)>29) towersize=0.175;
364  else{
365  if(abs(ieta)==22) towersize=0.1;
366  else if(abs(ieta)==23) towersize=0.113;
367  else if(abs(ieta)==24) towersize=0.129;
368  else if(abs(ieta)==25) towersize=0.16;
369  else if(abs(ieta)==26) towersize=0.168;
370  else if(abs(ieta)==27) towersize=0.15;
371  else if(abs(ieta)==28) towersize=0.218;
372  else if(abs(ieta)==29) towersize=0.132;
373  }
374  }
375 
376  deltaR /= towersize;
377  //if(deltaR<1/3.) nclosetowers++;
378 
379  if(deltaR<minDeltaR){
380  returnId = DetId(*did);
381  minDeltaR = deltaR;
382  }
383 
384  //if(abs(eta-hcalEta)<towersize/2. && abs(phi-hcalPhi)<towersize/2.) break;
385  }
386  }
387  //if(nclosetowers>1)std::cout<<"eta "<<eta<<" phi "<<phi<<" minDeltaR "<<minDeltaR<<" nclosetowers "<<nclosetowers<<std::endl;
388  return returnId;
389 
390 
391 }
392 
393 
395 
396  // get closest tower by delta R matching
397  // Now obsolute, instead use faster premade eta-phi grid
398 
399  double minDeltaR = 9999;
400 
401  std::vector<DetId> alldid = geo_->getValidDetIds();
402 
403  DetId returnId;
404 
405  //int nclosetowers=0;
406 
407  for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
408  if( (*did).det() == DetId::Hcal ){
409 
410  HcalDetId hid(*did);
411 
412  // which layer is irrelevant for an eta-phi map, no?
413 
414  if( hid.depth() != 1 ) continue;
415 
416  GlobalPoint pos =geo_->getGeometry(hid)->getPosition();
417 
418  double hcalEta = pos.eta();
419  double hcalPhi = pos.phi();
420 
421  //std::cout<<" transverse distance = "<<sqrt(pos.x()*pos.x()+pos.y()*pos.y())<<std::endl;
422 
423  //std::cout<<" ieta "<<(hid).ieta()<<" iphi "<<(hid).iphi()<<" hcalEta "<<hcalEta<<" hcalPhi "<<hcalPhi<<std::endl;
424 
425  double deltaR = reco::deltaR(eta,phi,hcalEta,hcalPhi);
426 
427  // need to factor in the size of the tower
428  double towersize = 0.087;
429 
430  int ieta = (hid).ieta();
431 
432  if(abs(ieta)>21){
433  if(abs(ieta)>29) towersize=0.175;
434  else{
435  if(abs(ieta)==22) towersize=0.1;
436  else if(abs(ieta)==23) towersize=0.113;
437  else if(abs(ieta)==24) towersize=0.129;
438  else if(abs(ieta)==25) towersize=0.16;
439  else if(abs(ieta)==26) towersize=0.168;
440  else if(abs(ieta)==27) towersize=0.15;
441  else if(abs(ieta)==28) towersize=0.218;
442  else if(abs(ieta)==29) towersize=0.132;
443  }
444  }
445 
446  deltaR /= towersize;
447  //if(deltaR<1/3.) nclosetowers++;
448 
449  if(deltaR<minDeltaR){
450  returnId = DetId(*did);
451  minDeltaR = deltaR;
452  }
453 
454  //if(abs(eta-hcalEta)<towersize/2. && abs(phi-hcalPhi)<towersize/2.) break;
455  }
456  }
457  //if(nclosetowers>1)std::cout<<"eta "<<eta<<" phi "<<phi<<" minDeltaR "<<minDeltaR<<" nclosetowers "<<nclosetowers<<std::endl;
458  return returnId;
459 
460 
461 }
462 
463 
464 
465 
466 // Taken from FastSimulation/CalorimeterProperties/src/HCALProperties.cc
467 // Note this returns an abs(ieta)
469  // binary search in the array of towers eta edges
470 
471  int size = 42;
472  if(!useHF_) size = 30;
473 
474  if(fabs(eta)>etaedge[size-1]) return -1;
475 
476  double x = fabs(eta);
477  int curr = size / 2;
478  int step = size / 4;
479  int iter;
480  int prevdir = 0;
481  int actudir = 0;
482 
483  for (iter = 0; iter < size ; iter++) {
484 
485  if( curr >= size || curr < 1 )
486  std::cout << " ParticleTowerProducer::eta2ieta - wrong current index = "
487  << curr << " !!!" << std::endl;
488 
489  if ((x <= etaedge[curr]) && (x > etaedge[curr-1])) break;
490  prevdir = actudir;
491  if(x > etaedge[curr]) {actudir = 1;}
492  else {actudir = -1;}
493  if(prevdir * actudir < 0) { if(step > 1) step /= 2;}
494  curr += actudir * step;
495  if(curr > size) curr = size;
496  else { if(curr < 1) {curr = 1;}}
497 
498  /*
499  std::cout << " HCALProperties::eta2ieta end of iter." << iter
500  << " curr, etaedge[curr-1], etaedge[curr] = "
501  << curr << " " << etaedge[curr-1] << " " << etaedge[curr] << std::endl;
502  */
503 
504  }
505 
506  /*
507  std::cout << " HCALProperties::eta2ieta for input x = " << x
508  << " found index = " << curr-1
509  << std::endl;
510  */
511 
512  return curr;
513 }
514 
515 int ParticleTowerProducer::phi2iphi(double phi, int ieta) const {
516 
517  if(phi<0) phi += 2.*PI;
518  else if(phi> 2.*PI) phi -= 2.*PI;
519 
520  int iphi = (int) TMath::Ceil(phi/2.0/PI*72.);
521  // take into account larger granularity in endcap (x2) and at the end of the HF (x4)
522  if(abs(ieta)>20){
523  if(abs(ieta)<40) iphi -= (iphi+1)%2;
524  else {
525  iphi -= (iphi+1)%4;
526  if(iphi==-1) iphi=71;
527  }
528  }
529 
530  return iphi;
531 
532 }
size
Write out results.
const double Pi
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
double eta() const final
momentum pseudorapidity
static const double etatow[]
int phi2iphi(double phi, int ieta) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
#define nullptr
ParticleTowerProducer(const edm::ParameterSet &)
void resetTowers(edm::Event &iEvent, const edm::EventSetup &iSetup)
int depth() const
get the tower depth
Definition: HcalDetId.cc:108
int iEvent
Definition: GenABIO.cc:230
std::map< DetId, double > towers_
double et() const final
transverse energy
double p4[4]
Definition: TauolaWrapper.h:92
void produce(edm::Event &, const edm::EventSetup &) override
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
static const double etacent[]
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
Definition: DetId.h:18
int eta2ieta(double eta) const
double sd
const T & get() const
Definition: EventSetup.h:58
edm::EDGetTokenT< reco::PFCandidateCollection > src_
et
define resolution functions of each parameter
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:97
T eta() const
Definition: PV3DBase.h:76
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
HLT enums.
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Particle.h:23
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:80
CaloGeometry const * geo_
DetId getNearestTower(const reco::PFCandidate &in) const
step
double phi() const final
momentum azimuthal angle
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:510