CMS 3D CMS Logo

PixelQuadrupletEDProducer.cc
Go to the documentation of this file.
10 
14 
16 #include "LayerQuadruplets.h"
17 
19 public:
21  ~PixelQuadrupletEDProducer() = default;
22 
23  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
24 
25  virtual void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
26 
27 private:
29 
31 
33 };
34 
36  tripletToken_(consumes<IntermediateHitTriplets>(iConfig.getParameter<edm::InputTag>("triplets"))),
37  generator_(iConfig, consumesCollector())
38 {
39  produces<RegionsSeedingHitSets>();
40 }
41 
44 
45  desc.add<edm::InputTag>("triplets", edm::InputTag("hitTripletEDProducer"));
47 
48  descriptions.add("pixelQuadrupletEDProducer", desc);
49 }
50 
53  iEvent.getByToken(tripletToken_, htriplets);
54  const auto& regionTriplets = *htriplets;
55 
56  const SeedingLayerSetsHits& seedingLayerHits = regionTriplets.seedingLayerHits();
57  if(seedingLayerHits.numberOfLayersInSet() < 4) {
58  throw cms::Exception("LogicError") << "PixelQuadrupletEDProducer expects SeedingLayerSetsHits::numberOfLayersInSet() to be >= 4, got " << seedingLayerHits.numberOfLayersInSet() << ". This is likely caused by a configuration error of this module, HitPairEDProducer, or SeedingLayersEDProducer.";
59  }
60 
61  auto seedingHitSets = std::make_unique<RegionsSeedingHitSets>();
62  if(regionTriplets.empty()) {
63  iEvent.put(std::move(seedingHitSets));
64  return;
65  }
66  seedingHitSets->reserve(regionTriplets.regionSize(), localRA_.upper());
67 
68  // match-making of triplet and quadruplet layers
69  std::vector<LayerQuadruplets::LayerSetAndLayers> quadlayers = LayerQuadruplets::layers(seedingLayerHits);
70 
71  LogDebug("PixelQuadrupletEDProducer") << "Creating quadruplets for " << regionTriplets.regionSize() << " regions, and " << quadlayers.size() << " triplet+4th layers from " << regionTriplets.tripletsSize() << " triplets";
72 
73  OrderedHitSeeds quadruplets;
74  quadruplets.reserve(localRA_.upper());
75 
76  for(const auto& regionLayerTriplets: regionTriplets) {
77  const TrackingRegion& region = regionLayerTriplets.region();
78  auto seedingHitSetsFiller = seedingHitSets->beginRegion(&region);
79 
80  LayerHitMapCache hitCache;
81  hitCache.extend(regionLayerTriplets.layerHitMapCache());
82 
83  LogTrace("PixelQuadrupletEDProducer") << " starting region, number of layer triplets " << regionLayerTriplets.layerTripletsSize();
84 
85  for(const auto& layerTriplet: regionLayerTriplets) {
86  LogTrace("PixelQuadrupletEDProducer") << " starting layer triplet " << layerTriplet.innerLayerIndex() << "," << layerTriplet.middleLayerIndex() << "," << layerTriplet.outerLayerIndex();
87  auto found = std::find_if(quadlayers.begin(), quadlayers.end(), [&](const LayerQuadruplets::LayerSetAndLayers& a) {
88  return a.first[0].index() == layerTriplet.innerLayerIndex() &&
89  a.first[1].index() == layerTriplet.middleLayerIndex() &&
90  a.first[2].index() == layerTriplet.outerLayerIndex();
91  });
92  if(found == quadlayers.end()) {
93  auto exp = cms::Exception("LogicError") << "Did not find the layer triplet from vector<triplet+fourth layers>. This is a sign of some internal inconsistency\n";
94  exp << "I was looking for layer triplet " << layerTriplet.innerLayerIndex() << "," << layerTriplet.middleLayerIndex() << "," << layerTriplet.outerLayerIndex()
95  << ". Quadruplets have the following triplets:\n";
96  for(const auto& a: quadlayers) {
97  exp << " " << a.first[0].index() << "," << a.first[1].index() << "," << a.first[2].index() << ": 4th layers";
98  for(const auto& b: a.second) {
99  exp << " " << b.index();
100  }
101  exp << "\n";
102  }
103  throw exp;
104  }
105  const auto& fourthLayers = found->second;
106 
107  generator_.hitQuadruplets(region, quadruplets, iEvent, iSetup, layerTriplet.begin(), layerTriplet.end(), fourthLayers, hitCache);
108 
109 #ifdef EDM_ML_DEBUG
110  LogTrace("PixelQuadrupletEDProducer") << " created " << quadruplets.size() << " quadruplets for layer triplet " << layerTriplet.innerLayerIndex() << "," << layerTriplet.middleLayerIndex() << "," << layerTriplet.outerLayerIndex() << " and 4th layers";
111  for(const auto& l: fourthLayers) {
112  LogTrace("PixelQuadrupletEDProducer") << " " << l.index();
113  }
114 #endif
115 
116  for(const auto& quad: quadruplets) {
117  seedingHitSetsFiller.emplace_back(quad[0], quad[1], quad[2], quad[3]);
118  }
119  quadruplets.clear();
120  }
121  }
122  localRA_.update(seedingHitSets->size());
123 
124  seedingHitSets->shrink_to_fit();
125  iEvent.put(std::move(seedingHitSets));
126 }
127 
#define LogDebug(id)
unsigned short numberOfLayersInSet() const
Get number of layers in each SeedingLayerSets.
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
std::pair< LayerSet, std::vector< Layer > > LayerSetAndLayers
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
static void fillDescriptions(edm::ParameterSetDescription &desc)
virtual unsigned int size() const
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int upper() const
int iEvent
Definition: GenABIO.cc:230
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
virtual void hitQuadruplets(const TrackingRegion &region, OrderedHitSeeds &result, const edm::Event &ev, const edm::EventSetup &es, const SeedingLayerSetsHits::SeedingLayerSet &tripletLayers, const std::vector< SeedingLayerSetsHits::SeedingLayer > &fourthLayers) override
edm::EDGetTokenT< IntermediateHitTriplets > tripletToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define LogTrace(id)
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
PixelQuadrupletEDProducer(const edm::ParameterSet &iConfig)
~PixelQuadrupletEDProducer()=default
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void extend(const LayerHitMapCache &other)
HLT enums.
double a
Definition: hdecay.h:121
PixelQuadrupletGenerator generator_
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:510
void update(unsigned int q)