CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
ParticleTowerProducer Class Reference

#include <RecoHI/ParticleTowerProducer/src/ParticleTowerProducer.cc>

Inheritance diagram for ParticleTowerProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 ParticleTowerProducer (const edm::ParameterSet &)
 
 ~ParticleTowerProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

virtual void beginJob ()
 
virtual void endJob ()
 
int eta2ieta (double eta) const
 
DetId getNearestTower (const reco::PFCandidate &in) const
 
DetId getNearestTower (double eta, double phi) const
 
int phi2iphi (double phi, int ieta) const
 
virtual void produce (edm::Event &, const edm::EventSetup &)
 
void resetTowers (edm::Event &iEvent, const edm::EventSetup &iSetup)
 

Private Attributes

double etaedge [42]
 
CaloGeometry const * geo_
 
double PI
 
TRandom * random_
 
edm::EDGetTokenT
< reco::PFCandidateCollection
src_
 
std::map< DetId, double > towers_
 
bool useHF_
 

Static Private Attributes

static const double etacent []
 
static const double etatow []
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 30 of file ParticleTowerProducer.h.

Constructor & Destructor Documentation

ParticleTowerProducer::ParticleTowerProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 65 of file ParticleTowerProducer.cc.

References edm::ParameterSet::getParameter(), Pi, PI, random_, src_, and useHF_.

65  :
66  geo_(0)
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 }
const double Pi
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::PFCandidateCollection > src_
CaloGeometry const * geo_
ParticleTowerProducer::~ParticleTowerProducer ( )

Definition at line 83 of file ParticleTowerProducer.cc.

84 {
85 
86  // do anything here that needs to be done at desctruction time
87  // (e.g. close files, deallocate resources etc.)
88 
89 }

Member Function Documentation

void ParticleTowerProducer::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDProducer.

Definition at line 211 of file ParticleTowerProducer.cc.

References etacent, etaedge, etatow, i, and useHF_.

212 {
213  // tower edges from fast sim, used starting at index 30 for the HF
214  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};
215 
216  // 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
217  const double etacent[40] = {
218  0.0435, // 1
219  0.1305, // 2
220  0.2175, // 3
221  0.3045, // 4
222  0.3915, // 5
223  0.4785, // 6
224  0.5655, // 7
225  0.6525, // 8
226  0.7395, // 9
227  0.8265, // 10
228  0.9135, // 11
229  1.0005, // 12
230  1.0875, // 13
231  1.1745, // 14
232  1.2615, // 15
233  1.3485, // 16
234  1.4355, // 17
235  1.5225, // 18
236  1.6095, // 19
237  1.6965, // 20
238  1.785 , // 21
239  1.88 , // 22
240  1.9865, // 23
241  2.1075, // 24
242  2.247 , // 25
243  2.411 , // 26
244  2.575 , // 27
245  2.759 , // 28
246  2.934 , // 29
247  2.90138,// 29
248  3.04448,// 30
249  3.21932,// 31
250  3.39462,// 32
251  3.56966,// 33
252  3.74485,// 34
253  3.91956,// 35
254  4.09497,// 36
255  4.27007,// 37
256  4.44417,// 38
257  4.62046 // 39
258  };
259 
260  // Use the real towers centrers for the barrel and endcap up to eta=3
261  etaedge[0] = 0.;
262  for(int i=1;i<30;i++){
263  etaedge[i] = (etacent[i-1]-etaedge[i-1])*2.0 + etaedge[i-1];
264  //std::cout<<" i "<<i<<" etaedge "<<etaedge[i]<<std::endl;
265  }
266  // beyond eta = 3 just use the fast sim values
267  if(useHF_){
268  for(int i=30;i<42;i++){
269  etaedge[i]=etatow[i];
270  //std::cout<<" i "<<i<<" etaedge "<<etaedge[i]<<std::endl;
271  }
272  }
273 }
int i
Definition: DBlmapReader.cc:9
static const double etatow[]
static const double etacent[]
void ParticleTowerProducer::endJob ( void  )
privatevirtual

Reimplemented from edm::EDProducer.

Definition at line 277 of file ParticleTowerProducer.cc.

277  {
278 }
int ParticleTowerProducer::eta2ieta ( double  eta) const
private

Definition at line 467 of file ParticleTowerProducer.cc.

References gather_cfg::cout, etaedge, findQualityFiles::size, relval_parameters_module::step, useHF_, and x.

Referenced by produce().

467  {
468  // binary search in the array of towers eta edges
469 
470  int size = 42;
471  if(!useHF_) size = 30;
472 
473  if(fabs(eta)>etaedge[size-1]) return -1;
474 
475  double x = fabs(eta);
476  int curr = size / 2;
477  int step = size / 4;
478  int iter;
479  int prevdir = 0;
480  int actudir = 0;
481 
482  for (iter = 0; iter < size ; iter++) {
483 
484  if( curr >= size || curr < 1 )
485  std::cout << " ParticleTowerProducer::eta2ieta - wrong current index = "
486  << curr << " !!!" << std::endl;
487 
488  if ((x <= etaedge[curr]) && (x > etaedge[curr-1])) break;
489  prevdir = actudir;
490  if(x > etaedge[curr]) {actudir = 1;}
491  else {actudir = -1;}
492  if(prevdir * actudir < 0) { if(step > 1) step /= 2;}
493  curr += actudir * step;
494  if(curr > size) curr = size;
495  else { if(curr < 1) {curr = 1;}}
496 
497  /*
498  std::cout << " HCALProperties::eta2ieta end of iter." << iter
499  << " curr, etaedge[curr-1], etaedge[curr] = "
500  << curr << " " << etaedge[curr-1] << " " << etaedge[curr] << std::endl;
501  */
502 
503  }
504 
505  /*
506  std::cout << " HCALProperties::eta2ieta for input x = " << x
507  << " found index = " << curr-1
508  << std::endl;
509  */
510 
511  return curr;
512 }
tuple cout
Definition: gather_cfg.py:121
tuple size
Write out results.
DetId ParticleTowerProducer::getNearestTower ( const reco::PFCandidate in) const
private

Definition at line 323 of file ParticleTowerProducer.cc.

References funct::abs(), deltaR(), reco::deltaR(), HcalDetId::depth(), eta, PV3DBase< T, PVType, FrameType >::eta(), reco::LeafCandidate::eta(), geo_, CaloGeometry::getGeometry(), CaloCellGeometry::getPosition(), CaloGeometry::getValidDetIds(), DetId::Hcal, HLT_25ns14e33_v1_cff::minDeltaR, phi, PV3DBase< T, PVType, FrameType >::phi(), and reco::LeafCandidate::phi().

323  {
324 
325  double eta = in.eta();
326  double phi = in.phi();
327 
328  double minDeltaR = 9999;
329 
330  std::vector<DetId> alldid = geo_->getValidDetIds();
331 
332  DetId returnId;
333 
334  //int nclosetowers=0;
335 
336  for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
337  if( (*did).det() == DetId::Hcal ){
338 
339  HcalDetId hid(*did);
340 
341  // which layer is irrelevant for an eta-phi map, no?
342 
343  if( hid.depth() != 1 ) continue;
344 
345  GlobalPoint pos =geo_->getGeometry(hid)->getPosition();
346 
347  double hcalEta = pos.eta();
348  double hcalPhi = pos.phi();
349 
350  //std::cout<<" transverse distance = "<<sqrt(pos.x()*pos.x()+pos.y()*pos.y())<<std::endl;
351 
352  //std::cout<<" ieta "<<(hid).ieta()<<" iphi "<<(hid).iphi()<<" hcalEta "<<hcalEta<<" hcalPhi "<<hcalPhi<<std::endl;
353 
354  double deltaR = reco::deltaR(eta,phi,hcalEta,hcalPhi);
355 
356  // need to factor in the size of the tower
357  double towersize = 0.087;
358 
359  int ieta = (hid).ieta();
360 
361  if(abs(ieta)>21){
362  if(abs(ieta)>29) towersize=0.175;
363  else{
364  if(abs(ieta)==22) towersize=0.1;
365  else if(abs(ieta)==23) towersize=0.113;
366  else if(abs(ieta)==24) towersize=0.129;
367  else if(abs(ieta)==25) towersize=0.16;
368  else if(abs(ieta)==26) towersize=0.168;
369  else if(abs(ieta)==27) towersize=0.15;
370  else if(abs(ieta)==28) towersize=0.218;
371  else if(abs(ieta)==29) towersize=0.132;
372  }
373  }
374 
375  deltaR /= towersize;
376  //if(deltaR<1/3.) nclosetowers++;
377 
378  if(deltaR<minDeltaR){
379  returnId = DetId(*did);
380  minDeltaR = deltaR;
381  }
382 
383  //if(abs(eta-hcalEta)<towersize/2. && abs(phi-hcalPhi)<towersize/2.) break;
384  }
385  }
386  //if(nclosetowers>1)std::cout<<"eta "<<eta<<" phi "<<phi<<" minDeltaR "<<minDeltaR<<" nclosetowers "<<nclosetowers<<std::endl;
387  return returnId;
388 
389 
390 }
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
virtual double eta() const
momentum pseudorapidity
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
Definition: DetId.h:18
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:90
T eta() const
Definition: PV3DBase.h:76
const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:76
CaloGeometry const * geo_
const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
virtual double phi() const
momentum azimuthal angle
DetId ParticleTowerProducer::getNearestTower ( double  eta,
double  phi 
) const
private

Definition at line 393 of file ParticleTowerProducer.cc.

References funct::abs(), deltaR(), reco::deltaR(), HcalDetId::depth(), PV3DBase< T, PVType, FrameType >::eta(), geo_, CaloGeometry::getGeometry(), CaloCellGeometry::getPosition(), CaloGeometry::getValidDetIds(), DetId::Hcal, HLT_25ns14e33_v1_cff::minDeltaR, and PV3DBase< T, PVType, FrameType >::phi().

393  {
394 
395  // get closest tower by delta R matching
396  // Now obsolute, instead use faster premade eta-phi grid
397 
398  double minDeltaR = 9999;
399 
400  std::vector<DetId> alldid = geo_->getValidDetIds();
401 
402  DetId returnId;
403 
404  //int nclosetowers=0;
405 
406  for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
407  if( (*did).det() == DetId::Hcal ){
408 
409  HcalDetId hid(*did);
410 
411  // which layer is irrelevant for an eta-phi map, no?
412 
413  if( hid.depth() != 1 ) continue;
414 
415  GlobalPoint pos =geo_->getGeometry(hid)->getPosition();
416 
417  double hcalEta = pos.eta();
418  double hcalPhi = pos.phi();
419 
420  //std::cout<<" transverse distance = "<<sqrt(pos.x()*pos.x()+pos.y()*pos.y())<<std::endl;
421 
422  //std::cout<<" ieta "<<(hid).ieta()<<" iphi "<<(hid).iphi()<<" hcalEta "<<hcalEta<<" hcalPhi "<<hcalPhi<<std::endl;
423 
424  double deltaR = reco::deltaR(eta,phi,hcalEta,hcalPhi);
425 
426  // need to factor in the size of the tower
427  double towersize = 0.087;
428 
429  int ieta = (hid).ieta();
430 
431  if(abs(ieta)>21){
432  if(abs(ieta)>29) towersize=0.175;
433  else{
434  if(abs(ieta)==22) towersize=0.1;
435  else if(abs(ieta)==23) towersize=0.113;
436  else if(abs(ieta)==24) towersize=0.129;
437  else if(abs(ieta)==25) towersize=0.16;
438  else if(abs(ieta)==26) towersize=0.168;
439  else if(abs(ieta)==27) towersize=0.15;
440  else if(abs(ieta)==28) towersize=0.218;
441  else if(abs(ieta)==29) towersize=0.132;
442  }
443  }
444 
445  deltaR /= towersize;
446  //if(deltaR<1/3.) nclosetowers++;
447 
448  if(deltaR<minDeltaR){
449  returnId = DetId(*did);
450  minDeltaR = deltaR;
451  }
452 
453  //if(abs(eta-hcalEta)<towersize/2. && abs(phi-hcalPhi)<towersize/2.) break;
454  }
455  }
456  //if(nclosetowers>1)std::cout<<"eta "<<eta<<" phi "<<phi<<" minDeltaR "<<minDeltaR<<" nclosetowers "<<nclosetowers<<std::endl;
457  return returnId;
458 
459 
460 }
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
Definition: DetId.h:18
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:90
T eta() const
Definition: PV3DBase.h:76
const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:76
CaloGeometry const * geo_
const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
int ParticleTowerProducer::phi2iphi ( double  phi,
int  ieta 
) const
private

Definition at line 514 of file ParticleTowerProducer.cc.

References funct::abs(), and PI.

Referenced by produce().

514  {
515 
516  if(phi<0) phi += 2.*PI;
517  else if(phi> 2.*PI) phi -= 2.*PI;
518 
519  int iphi = (int) TMath::Ceil(phi/2.0/PI*72.);
520  // take into account larger granularity in endcap (x2) and at the end of the HF (x4)
521  if(abs(ieta)>20){
522  if(abs(ieta)<40) iphi -= (iphi+1)%2;
523  else {
524  iphi -= (iphi+1)%4;
525  if(iphi==-1) iphi=71;
526  }
527  }
528 
529  return iphi;
530 
531 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void ParticleTowerProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
privatevirtual

Implements edm::EDProducer.

Definition at line 98 of file ParticleTowerProducer.cc.

References funct::abs(), reco::LeafCandidate::et(), eta, PV3DBase< T, PVType, FrameType >::eta(), reco::LeafCandidate::eta(), eta2ieta(), geo_, edm::EventSetup::get(), edm::Event::getByToken(), CaloGeometry::getGeometry(), CaloCellGeometry::getPosition(), HcalBarrel, HcalEndcap, HcalForward, p4, PV3DBase< T, PVType, FrameType >::phi(), reco::LeafCandidate::phi(), phi2iphi(), parseEventContent::prod, edm::ESHandle< class >::product(), edm::Event::put(), resetTowers(), sd, src_, towers_, and useHF_.

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  std::auto_ptr<CaloTowerCollection> prod(new CaloTowerCollection());
164 
165  for ( std::map< DetId, double >::const_iterator iter = towers_.begin();
166  iter != towers_.end(); ++iter ){
167 
168  CaloTowerDetId newTowerId(iter->first.rawId());
169  double et = iter->second;
170 
171  if(et>0){
172 
173  GlobalPoint pos =geo_->getGeometry(newTowerId)->getPosition();
174 
175  if(!useHF_ && fabs(pos.eta()) > 3. ) continue;
176 
177  // currently sets et = pt, mass to zero
178  // pt, eta , phi, mass
179  reco::Particle::PolarLorentzVector p4(et,pos.eta(),pos.phi(),0.);
180 
181  CaloTower newTower(newTowerId,et,0,0,0,0,p4,pos,pos);
182  prod->push_back(newTower);
183  }
184  }
185 
186 
187  //For reference, Calo Tower Constructors
188 
189  /*
190  CaloTower(const CaloTowerDetId& id,
191  double emE, double hadE, double outerE,
192  int ecal_tp, int hcal_tp,
193  const PolarLorentzVector p4,
194  GlobalPoint emPosition, GlobalPoint hadPosition);
195 
196  CaloTower(const CaloTowerDetId& id,
197  double emE, double hadE, double outerE,
198  int ecal_tp, int hcal_tp,
199  const LorentzVector p4,
200  GlobalPoint emPosition, GlobalPoint hadPosition);
201  */
202 
203 
204  iEvent.put(prod);
205 
206 
207 }
int phi2iphi(double phi, int ieta) const
virtual double et() const
transverse energy
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
virtual double eta() const
momentum pseudorapidity
void resetTowers(edm::Event &iEvent, const edm::EventSetup &iSetup)
std::map< DetId, double > towers_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
double p4[4]
Definition: TauolaWrapper.h:92
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int eta2ieta(double eta) const
double sd
const T & get() const
Definition: EventSetup.h:56
edm::EDGetTokenT< reco::PFCandidateCollection > src_
T const * product() const
Definition: ESHandle.h:86
T eta() const
Definition: PV3DBase.h:76
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
edm::SortedCollection< CaloTower > CaloTowerCollection
Definition: CaloTowerFwd.h:15
const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:76
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Particle.h:23
CaloGeometry const * geo_
const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
virtual double phi() const
momentum azimuthal angle
void ParticleTowerProducer::resetTowers ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
private

Definition at line 281 of file ParticleTowerProducer.cc.

References HcalDetId::depth(), PV3DBase< T, PVType, FrameType >::eta(), geo_, CaloGeometry::getGeometry(), CaloCellGeometry::getPosition(), CaloGeometry::getValidDetIds(), DetId::Hcal, towers_, and useHF_.

Referenced by produce().

282 {
283 
284  std::vector<DetId> alldid = geo_->getValidDetIds();
285 
286  for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
287  if( (*did).det() == DetId::Hcal ){
288  HcalDetId hid = HcalDetId(*did);
289  if( hid.depth() == 1 ) {
290 
291  if(!useHF_){
292  GlobalPoint pos =geo_->getGeometry(hid)->getPosition();
293  //if((hid).iphi()==1)std::cout<<" ieta "<<(hid).ieta()<<" eta "<<pos.eta()<<" iphi "<<(hid).iphi()<<" phi "<<pos.phi()<<std::endl;
294  if(fabs(pos.eta())>3.) continue;
295  }
296  towers_[(*did)] = 0.;
297  }
298 
299  }
300  }
301 
302 
303 
304 
305  //print out tower eta boundaries
306  /*
307  int current_ieta=0;
308  float testphi =0.;
309  for(double testeta=1.7; testeta<=5.4;testeta+=0.00001){
310  HcalDetId hid = getNearestTower(testeta,testphi);
311  if((hid).ieta()!=current_ieta) {
312  std::cout<<" new ieta "<<current_ieta<<" testeta "<<testeta<<std::endl;
313  current_ieta = (hid).ieta();
314  }
315  }
316  */
317 
318 }
int depth() const
get the tower depth
Definition: HcalDetId.h:55
std::map< DetId, double > towers_
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:90
T eta() const
Definition: PV3DBase.h:76
const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:76
CaloGeometry const * geo_
const GlobalPoint & getPosition() const
Returns the position of reference for this cell.

Member Data Documentation

const double ParticleTowerProducer::etacent[]
staticprivate

Definition at line 61 of file ParticleTowerProducer.h.

Referenced by beginJob().

double ParticleTowerProducer::etaedge[42]
private

Definition at line 62 of file ParticleTowerProducer.h.

Referenced by beginJob(), and eta2ieta().

const double ParticleTowerProducer::etatow[]
staticprivate

Definition at line 60 of file ParticleTowerProducer.h.

Referenced by beginJob().

CaloGeometry const* ParticleTowerProducer::geo_
private

Definition at line 57 of file ParticleTowerProducer.h.

Referenced by getNearestTower(), produce(), and resetTowers().

double ParticleTowerProducer::PI
private

Definition at line 54 of file ParticleTowerProducer.h.

Referenced by ParticleTowerProducer(), and phi2iphi().

TRandom* ParticleTowerProducer::random_
private

Definition at line 55 of file ParticleTowerProducer.h.

Referenced by ParticleTowerProducer().

edm::EDGetTokenT<reco::PFCandidateCollection> ParticleTowerProducer::src_
private

Definition at line 48 of file ParticleTowerProducer.h.

Referenced by ParticleTowerProducer(), and produce().

std::map<DetId,double> ParticleTowerProducer::towers_
private

Definition at line 51 of file ParticleTowerProducer.h.

Referenced by produce(), and resetTowers().

bool ParticleTowerProducer::useHF_
private