CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MultiHitFromChi2EDProducer.cc
Go to the documentation of this file.
10 
16 
18 public:
20  ~MultiHitFromChi2EDProducer() override = default;
21 
22  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
23 
24  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
25 
26 private:
28 
30 
32 };
33 
35  : doubletToken_(consumes<IntermediateHitDoublets>(iConfig.getParameter<edm::InputTag>("doublets"))),
36  generator_(iConfig, consumesCollector()) {
37  produces<RegionsSeedingHitSets>();
38  produces<edm::OwnVector<BaseTrackerRecHit> >();
39 }
40 
43 
44  desc.add<edm::InputTag>("doublets", edm::InputTag("hitPairEDProducer"));
45 
47 
49  descriptions.add(label, desc);
50 }
51 
54  iEvent.getByToken(doubletToken_, hdoublets);
55  const auto& regionDoublets = *hdoublets;
56 
57  const SeedingLayerSetsHits& seedingLayerHits = regionDoublets.seedingLayerHits();
58  if (seedingLayerHits.numberOfLayersInSet() < 3) {
59  throw cms::Exception("LogicError")
60  << "MultiHitFromChi2EDProducer expects SeedingLayerSetsHits::numberOfLayersInSet() to be >= 3, got "
61  << seedingLayerHits.numberOfLayersInSet()
62  << ". This is likely caused by a configuration error of this module, HitPairEDProducer, or "
63  "SeedingLayersEDProducer.";
64  }
65 
66  auto seedingHitSets = std::make_unique<RegionsSeedingHitSets>();
67  if (regionDoublets.empty()) {
68  iEvent.put(std::move(seedingHitSets));
69  return;
70  }
71  seedingHitSets->reserve(regionDoublets.regionSize(), localRA_.upper());
72  generator_.initES(iSetup);
73 
74  // match-making of pair and triplet layers
75  std::vector<LayerTriplets::LayerSetAndLayers> trilayers = LayerTriplets::layers(seedingLayerHits);
76 
77  OrderedMultiHits multihits;
78  multihits.reserve(localRA_.upper());
79  std::vector<std::unique_ptr<BaseTrackerRecHit> > refittedHitStorage;
80  refittedHitStorage.reserve(localRA_.upper() * 2);
81 
82  LogDebug("MultiHitFromChi2EDProducer") << "Creating multihits for " << regionDoublets.regionSize() << " regions, and "
83  << trilayers.size() << " pair+3rd layers from "
84  << regionDoublets.layerPairsSize() << " layer pairs";
85 
86  LayerHitMapCache hitCache;
87  for (const auto& regionLayerPairs : regionDoublets) {
88  const TrackingRegion& region = regionLayerPairs.region();
89 
90  auto seedingHitSetsFiller = RegionsSeedingHitSets::dummyFiller();
91  seedingHitSetsFiller = seedingHitSets->beginRegion(&region);
92 
93  hitCache.clear();
94  hitCache.extend(regionLayerPairs.layerHitMapCache());
95 
96  LogTrace("MultiHitFromChi2EDProducer") << " starting region";
97 
98  for (const auto& layerPair : regionLayerPairs) {
99  LogTrace("MultiHitFromChi2EDProducer")
100  << " starting layer pair " << layerPair.innerLayerIndex() << "," << layerPair.outerLayerIndex();
101 
102  auto found = std::find_if(trilayers.begin(), trilayers.end(), [&](const LayerTriplets::LayerSetAndLayers& a) {
103  return a.first[0].index() == layerPair.innerLayerIndex() && a.first[1].index() == layerPair.outerLayerIndex();
104  });
105  if (found == trilayers.end()) {
106  auto exp = cms::Exception("LogicError") << "Did not find the layer pair from vector<pair+third layers>. This "
107  "is a sign of some internal inconsistency\n";
108  exp << "I was looking for layer pair " << layerPair.innerLayerIndex() << "," << layerPair.outerLayerIndex()
109  << ". Triplets have the following pairs:\n";
110  for (const auto& a : trilayers) {
111  exp << " " << a.first[0].index() << "," << a.first[1].index() << ": 3rd layers";
112  for (const auto& b : a.second) {
113  exp << " " << b.index();
114  }
115  exp << "\n";
116  }
117  throw exp;
118  }
119  const auto& thirdLayers = found->second;
120 
121  generator_.hitSets(region, multihits, layerPair.doublets(), thirdLayers, hitCache, refittedHitStorage);
122 
123 #ifdef EDM_ML_DEBUG
124  LogTrace("MultiHitFromChi2EDProducer")
125  << " created " << multihits.size() << " multihits for layer pair " << layerPair.innerLayerIndex() << ","
126  << layerPair.outerLayerIndex() << " and 3rd layers";
127  for (const auto& l : thirdLayers) {
128  LogTrace("MultiHitFromChi2EDProducer") << " " << l.index();
129  }
130 #endif
131 
132  for (const SeedingHitSet& hitSet : multihits) {
133  seedingHitSetsFiller.emplace_back(hitSet);
134  }
135  multihits.clear();
136  }
137  }
138  localRA_.update(seedingHitSets->size());
139 
140  auto storage = std::make_unique<edm::OwnVector<BaseTrackerRecHit> >();
141  storage->reserve(refittedHitStorage.size());
142  for (auto& ptr : refittedHitStorage)
143  storage->push_back(ptr.release());
144 
145  seedingHitSets->shrink_to_fit();
146  storage->shrink_to_fit();
147  iEvent.put(std::move(seedingHitSets));
148  iEvent.put(std::move(storage));
149 }
150 
tuple seedingHitSets
unsigned short numberOfLayersInSet() const
Get number of layers in each SeedingLayerSets.
unsigned int size() const override
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
#define LogTrace(id)
char const * label
int upper() const
~MultiHitFromChi2EDProducer() override=default
int iEvent
Definition: GenABIO.cc:224
static RegionFiller dummyFiller()
def move
Definition: eostools.py:511
void initES(const edm::EventSetup &es) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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:118
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static const char * fillDescriptionsLabel()
void extend(const LayerHitMapCache &other)
double a
Definition: hdecay.h:119
static void fillDescriptions(edm::ParameterSetDescription &desc)
edm::EDGetTokenT< IntermediateHitDoublets > doubletToken_
void update(unsigned int q)
#define LogDebug(id)
MultiHitFromChi2EDProducer(const edm::ParameterSet &iConfig)