CMS 3D CMS Logo

CombinedHitPairGeneratorForPhotonConversion.cc
Go to the documentation of this file.
6 
8 namespace {
9  edm::RunningAverage localRA;
10 }
11 
12 
14  theSeedingLayerToken(iC.consumes<SeedingLayerSetsHits>(cfg.getParameter<edm::InputTag>("SeedingLayers")))
15 {
16  theMaxElement = cfg.getParameter<unsigned int>("maxElement");
17  maxHitPairsPerTrackAndGenerator = cfg.getParameter<unsigned int>("maxHitPairsPerTrackAndGenerator");
18  theGenerator.reset(new HitPairGeneratorFromLayerPairForPhotonConversion(0, 1, &theLayerCache, 0, maxHitPairsPerTrackAndGenerator));
19 }
20 
21 
22 
24  const ConversionRegion& convRegion,
25  const TrackingRegion& region, const edm::Event & ev, const edm::EventSetup& es)
26 {
27  if (thePairs.capacity()==0) thePairs.reserve(localRA.upper());
28  thePairs.clear();
29  hitPairs(convRegion, region, thePairs, ev, es);
30  return thePairs;
31 }
32 
33 
35  const ConversionRegion& convRegion,
36  const TrackingRegion& region, OrderedHitPairs & result,
37  const edm::Event& ev, const edm::EventSetup& es)
38 {
40  ev.getByToken(theSeedingLayerToken, hlayers);
41  const SeedingLayerSetsHits& layers = *hlayers;
42  assert(layers.numberOfLayersInSet() == 2);
43 
44  for(SeedingLayerSetsHits::LayerSetIndex i=0; i<layers.size(); ++i) {
45  theGenerator->hitPairs( convRegion, region, result, layers[i], ev, es);
46  }
47 
48 }
49 
50 
53  localRA.update(thePairs.size()); thePairs.clear(); thePairs.shrink_to_fit();
54 }
55 
T getParameter(std::string const &) const
unsigned short numberOfLayersInSet() const
Get number of layers in each SeedingLayerSets.
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
unsigned short LayerSetIndex
bool ev
const OrderedHitPairs & run(const ConversionRegion &convRegion, const TrackingRegion &region, const edm::Event &ev, const edm::EventSetup &es)
unsigned int size() const override
CombinedHitPairGeneratorForPhotonConversion(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
void hitPairs(const ConversionRegion &convRegion, const TrackingRegion &reg, OrderedHitPairs &result, const edm::Event &ev, const edm::EventSetup &es)
HLT enums.
std::unique_ptr< HitPairGeneratorFromLayerPairForPhotonConversion > theGenerator
unsigned short size() const
Get the number of SeedingLayerSets.