CMS 3D CMS Logo

HGC3DClusterGenMatchSelector.cc
Go to the documentation of this file.
6 
10 
11 namespace l1t {
13  public:
16 
17  private:
20  double dR_;
21  void produce(edm::Event &, const edm::EventSetup &) override;
22 
23  }; // class
24 } // namespace l1t
25 
27  : src_(consumes<l1t::HGCalMulticlusterBxCollection>(iConfig.getParameter<edm::InputTag>("src"))),
28  genParticleSrc_(consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("genSrc"))),
29  dR_(iConfig.getParameter<double>("dR")) {
30  produces<l1t::HGCalMulticlusterBxCollection>();
31 }
32 
34  auto out = std::make_unique<l1t::HGCalMulticlusterBxCollection>();
35 
37  iEvent.getByToken(src_, multiclusters);
38 
40  iEvent.getByToken(genParticleSrc_, genParticles);
41 
42  for (int bx = multiclusters->getFirstBX(); bx <= multiclusters->getLastBX(); ++bx) {
43  for (auto it = multiclusters->begin(bx), ed = multiclusters->end(bx); it != ed; ++it) {
44  const auto &multicluster = *it;
45  for (const auto &particle : *genParticles) {
46  if (particle.status() != 1)
47  continue;
48  if (deltaR(multicluster, particle) < dR_) {
49  out->push_back(bx, multicluster);
50  break; // don't save duplicate multiclusters!
51  }
52  }
53  }
54  }
55 
56  iEvent.put(std::move(out));
57 }
genParticles2HepMC_cfi.genParticles
genParticles
Definition: genParticles2HepMC_cfi.py:4
l1t::HGC3DClusterGenMatchSelector
Definition: HGC3DClusterGenMatchSelector.cc:12
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
reco::GenParticleCollection
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
Definition: GenParticleFwd.h:13
l1GtPatternGenerator_cfi.bx
bx
Definition: l1GtPatternGenerator_cfi.py:18
EDProducer.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:46
edm::Handle
Definition: AssociativeIterator.h:50
l1t::HGC3DClusterGenMatchSelector::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: HGC3DClusterGenMatchSelector.cc:33
GenParticle.h
BXVector
Definition: BXVector.h:15
BXVector::getFirstBX
int getFirstBX() const
MakerMacros.h
HGCalMulticluster.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
l1t::HGC3DClusterGenMatchSelector::genParticleSrc_
edm::EDGetToken genParticleSrc_
Definition: HGC3DClusterGenMatchSelector.cc:19
BXVector::begin
const_iterator begin(int bx) const
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
l1t::HGC3DClusterGenMatchSelector::HGC3DClusterGenMatchSelector
HGC3DClusterGenMatchSelector(const edm::ParameterSet &)
Definition: HGC3DClusterGenMatchSelector.cc:26
BXVector::end
const_iterator end(int bx) const
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
deltaR.h
l1t
delete x;
Definition: CaloConfig.h:22
l1t::HGC3DClusterGenMatchSelector::src_
edm::EDGetTokenT< l1t::HGCalMulticlusterBxCollection > src_
Definition: HGC3DClusterGenMatchSelector.cc:18
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::stream::EDProducer
Definition: EDProducer.h:36
edm::EventSetup
Definition: EventSetup.h:58
edm::EDGetToken
Definition: EDGetToken.h:35
l1t::HGC3DClusterGenMatchSelector::~HGC3DClusterGenMatchSelector
~HGC3DClusterGenMatchSelector() override
Definition: HGC3DClusterGenMatchSelector.cc:15
eostools.move
def move(src, dest)
Definition: eostools.py:511
Frameworkfwd.h
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
ParameterSet.h
edm::Event
Definition: Event.h:73
BXVector::getLastBX
int getLastBX() const
l1t::HGC3DClusterGenMatchSelector::dR_
double dR_
Definition: HGC3DClusterGenMatchSelector.cc:20