CMS 3D CMS Logo

MultiHitFromChi2EDProducer.cc
Go to the documentation of this file.
11 
16 
17 #include <memory>
18 #include <vector>
19 
21 public:
23  ~MultiHitFromChi2EDProducer() override = default;
24 
25  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
26 
27  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
28 
29 private:
32 
34 
36 };
37 
39  : doubletToken_(consumes<IntermediateHitDoublets>(iConfig.getParameter<edm::InputTag>("doublets"))),
40  hitsPutToken_{produces()},
41  generator_(iConfig, consumesCollector()) {
42  produces<RegionsSeedingHitSets>();
43 }
44 
47 
48  desc.add<edm::InputTag>("doublets", edm::InputTag("hitPairEDProducer"));
49 
51 
53  descriptions.add(label, desc);
54 }
55 
58  iEvent.getByToken(doubletToken_, hdoublets);
59  const auto& regionDoublets = *hdoublets;
60 
61  const SeedingLayerSetsHits& seedingLayerHits = regionDoublets.seedingLayerHits();
62  if (seedingLayerHits.numberOfLayersInSet() < 3) {
63  throw cms::Exception("LogicError")
64  << "MultiHitFromChi2EDProducer expects SeedingLayerSetsHits::numberOfLayersInSet() to be >= 3, got "
65  << seedingLayerHits.numberOfLayersInSet()
66  << ". This is likely caused by a configuration error of this module, HitPairEDProducer, or "
67  "SeedingLayersEDProducer.";
68  }
69 
70  auto seedingHitSets = std::make_unique<RegionsSeedingHitSets>();
71  if (regionDoublets.empty()) {
73  return;
74  }
75  seedingHitSets->reserve(regionDoublets.regionSize(), localRA_.upper());
76  generator_.initES(iSetup);
77 
78  // match-making of pair and triplet layers
79  std::vector<LayerTriplets::LayerSetAndLayers> trilayers = LayerTriplets::layers(seedingLayerHits);
80 
81  OrderedMultiHits multihits;
82  multihits.reserve(localRA_.upper());
83  std::vector<std::unique_ptr<BaseTrackerRecHit>> refittedHitStorage;
84  refittedHitStorage.reserve(localRA_.upper() * 2);
85 
86  LogDebug("MultiHitFromChi2EDProducer") << "Creating multihits for " << regionDoublets.regionSize() << " regions, and "
87  << trilayers.size() << " pair+3rd layers from "
88  << regionDoublets.layerPairsSize() << " layer pairs";
89 
90  LayerHitMapCache hitCache;
91  for (const auto& regionLayerPairs : regionDoublets) {
92  const TrackingRegion& region = regionLayerPairs.region();
93 
94  auto seedingHitSetsFiller = RegionsSeedingHitSets::dummyFiller();
95  seedingHitSetsFiller = seedingHitSets->beginRegion(&region);
96 
97  hitCache.clear();
98  hitCache.extend(regionLayerPairs.layerHitMapCache());
99 
100  LogTrace("MultiHitFromChi2EDProducer") << " starting region";
101 
102  for (const auto& layerPair : regionLayerPairs) {
103  LogTrace("MultiHitFromChi2EDProducer")
104  << " starting layer pair " << layerPair.innerLayerIndex() << "," << layerPair.outerLayerIndex();
105 
106  auto found = std::find_if(trilayers.begin(), trilayers.end(), [&](const LayerTriplets::LayerSetAndLayers& a) {
107  return a.first[0].index() == layerPair.innerLayerIndex() && a.first[1].index() == layerPair.outerLayerIndex();
108  });
109  if (found == trilayers.end()) {
110  auto exp = cms::Exception("LogicError") << "Did not find the layer pair from vector<pair+third layers>. This "
111  "is a sign of some internal inconsistency\n";
112  exp << "I was looking for layer pair " << layerPair.innerLayerIndex() << "," << layerPair.outerLayerIndex()
113  << ". Triplets have the following pairs:\n";
114  for (const auto& a : trilayers) {
115  exp << " " << a.first[0].index() << "," << a.first[1].index() << ": 3rd layers";
116  for (const auto& b : a.second) {
117  exp << " " << b.index();
118  }
119  exp << "\n";
120  }
121  throw exp;
122  }
123  const auto& thirdLayers = found->second;
124 
125  generator_.hitSets(region, multihits, layerPair.doublets(), thirdLayers, hitCache, refittedHitStorage);
126 
127 #ifdef EDM_ML_DEBUG
128  LogTrace("MultiHitFromChi2EDProducer")
129  << " created " << multihits.size() << " multihits for layer pair " << layerPair.innerLayerIndex() << ","
130  << layerPair.outerLayerIndex() << " and 3rd layers";
131  for (const auto& l : thirdLayers) {
132  LogTrace("MultiHitFromChi2EDProducer") << " " << l.index();
133  }
134 #endif
135 
136  for (const SeedingHitSet& hitSet : multihits) {
137  seedingHitSetsFiller.emplace_back(hitSet);
138  }
139  multihits.clear();
140  }
141  }
142  localRA_.update(seedingHitSets->size());
143 
144  seedingHitSets->shrink_to_fit();
145  refittedHitStorage.shrink_to_fit();
147  iEvent.emplace(hitsPutToken_, std::move(refittedHitStorage));
148 }
149 
unsigned int size() const override
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
edm::EDPutTokenT< std::vector< std::unique_ptr< BaseTrackerRecHit > > > hitsPutToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define LogTrace(id)
char const * label
~MultiHitFromChi2EDProducer() override=default
int iEvent
Definition: GenABIO.cc:224
static RegionFiller dummyFiller()
void initES(const edm::EventSetup &es) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
unsigned short numberOfLayersInSet() const
Get number of layers in each SeedingLayerSets.
std::pair< LayerSet, std::vector< Layer > > LayerSetAndLayers
Definition: LayerTriplets.h:16
MultiHitGeneratorFromChi2 generator_
void hitSets(const TrackingRegion &region, OrderedMultiHits &trs, const edm::Event &ev, const edm::EventSetup &es, SeedingLayerSetsHits::SeedingLayerSet pairLayers, std::vector< SeedingLayerSetsHits::SeedingLayer > thirdLayers) override
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static const char * fillDescriptionsLabel()
void extend(const LayerHitMapCache &other)
HLT enums.
double a
Definition: hdecay.h:121
static void fillDescriptions(edm::ParameterSetDescription &desc)
edm::EDGetTokenT< IntermediateHitDoublets > doubletToken_
def move(src, dest)
Definition: eostools.py:511
void update(unsigned int q)
#define LogDebug(id)
MultiHitFromChi2EDProducer(const edm::ParameterSet &iConfig)