CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaHLTCaloTowerProducer.cc
Go to the documentation of this file.
1 // makes CaloTowerCandidates from CaloTowers
2 // original author: M. Sani (UCSD)
3 
4 #include <cmath>
11 
14 
17 
18 // Math
19 #include "Math/GenVector/VectorUtil.h"
20 #include <cmath>
21 
22 using namespace edm;
23 using namespace reco;
24 using namespace std;
25 using namespace l1extra ;
26 
27 EgammaHLTCaloTowerProducer::EgammaHLTCaloTowerProducer( const ParameterSet & p ) : towers_ (consumes<CaloTowerCollection>(p.getParameter<InputTag> ("towerCollection"))),
28  cone_ (p.getParameter<double> ("useTowersInCone")),
29  l1isoseeds_ (consumes<l1extra::L1EmParticleCollection>(p.getParameter< edm::InputTag > ("L1IsoCand"))),
30  l1nonisoseeds_ (consumes<l1extra::L1EmParticleCollection>(p.getParameter< edm::InputTag > ("L1NonIsoCand"))),
31  EtThreshold_ (p.getParameter<double> ("EtMin")),
32  EThreshold_ (p.getParameter<double> ("EMin")) {
33 
34  produces<CaloTowerCollection>();
35 }
36 
38 
40 
41  desc.add<edm::InputTag>(("towerCollection"), edm::InputTag("hltRecoEcalCandidate"));
42  desc.add<edm::InputTag>(("L1IsoCand"), edm::InputTag("hltTowerMakerForAll"));
43  desc.add<edm::InputTag>(("L1NonIsoCand"), edm::InputTag("fixedGridRhoFastjetAllCalo"));
44  desc.add<double>(("useTowersInCone"), 0.8);
45  desc.add<double>(("EtMin"), 1.0);
46  desc.add<double>(("EMin"), 1.0);
47  descriptions.add(("hltCaloTowerForEgamma"), desc);
48 }
49 
50 
52 
54  evt.getByToken(towers_, caloTowers);
55 
57  evt.getByToken(l1isoseeds_, emIsolColl);
59  evt.getByToken(l1nonisoseeds_, emNonIsolColl);
60  std::auto_ptr<CaloTowerCollection> cands(new CaloTowerCollection);
61  cands->reserve(caloTowers->size());
62 
63  for (unsigned idx = 0; idx < caloTowers->size(); idx++) {
64  const CaloTower* cal = &((*caloTowers) [idx]);
65  if (cal->et() >= EtThreshold_ && cal->energy() >= EThreshold_) {
66  bool fill = false;
67  math::PtEtaPhiELorentzVector p(cal->et(), cal->eta(), cal->phi(), cal->energy());
68  for (l1extra::L1EmParticleCollection::const_iterator emItr = emIsolColl->begin(); emItr != emIsolColl->end() ;++emItr) {
69  double delta = ROOT::Math::VectorUtil::DeltaR((*emItr).p4().Vect(), p);
70  if(delta < cone_) {
71  cands->push_back(*cal);
72  fill = true;
73  break;
74  }
75  }
76 
77  if (!fill) {
78  for(l1extra::L1EmParticleCollection::const_iterator emItr = emNonIsolColl->begin(); emItr != emNonIsolColl->end() ;++emItr) {
79  double delta = ROOT::Math::VectorUtil::DeltaR((*emItr).p4().Vect(), p);
80  if(delta < cone_) {
81  cands->push_back(*cal);
82  break;
83  }
84  }
85  }
86  }
87  }
88 
89  evt.put( cands );
90 }
dbl * delta
Definition: mlp_gen.cc:36
edm::EDGetTokenT< l1extra::L1EmParticleCollection > l1nonisoseeds_
edm::EDGetTokenT< CaloTowerCollection > towers_
string fill
Definition: lumiContext.py:319
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual float phi() const
momentum azimuthal angle
EgammaHLTCaloTowerProducer(const edm::ParameterSet &)
virtual double energy() const
energy
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
virtual float eta() const
momentum pseudorapidity
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< l1extra::L1EmParticleCollection > l1isoseeds_
PtEtaPhiELorentzVectorD PtEtaPhiELorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:27
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void produce(edm::Event &e, const edm::EventSetup &) override
double et(double vtxZ) const
Definition: CaloTower.h:116
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)