File indexing completed on 2024-04-06 12:28:53
0001 #include "CombinedMultiHitGenerator.h"
0002
0003 #include "RecoTracker/TkHitPairs/interface/HitPairGeneratorFromLayerPair.h"
0004 #include "RecoTracker/TkSeedGenerator/interface/MultiHitGeneratorFromPairAndLayers.h"
0005 #include "RecoTracker/TkSeedGenerator/interface/MultiHitGeneratorFromPairAndLayersFactory.h"
0006 #include "RecoTracker/PixelSeeding/interface/LayerTriplets.h"
0007 #include "FWCore/Framework/interface/ConsumesCollector.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "DataFormats/Common/interface/Handle.h"
0010
0011 CombinedMultiHitGenerator::CombinedMultiHitGenerator(const edm::ParameterSet& cfg, edm::ConsumesCollector& iC)
0012 : theSeedingLayerToken(iC.consumes<SeedingLayerSetsHits>(cfg.getParameter<edm::InputTag>("SeedingLayers"))) {
0013 edm::ParameterSet generatorPSet = cfg.getParameter<edm::ParameterSet>("GeneratorPSet");
0014 std::string generatorName = generatorPSet.getParameter<std::string>("ComponentName");
0015 theGenerator = MultiHitGeneratorFromPairAndLayersFactory::get()->create(generatorName, generatorPSet, iC);
0016 theGenerator->init(std::make_unique<HitPairGeneratorFromLayerPair>(iC, 0, 1, &theLayerCache), &theLayerCache);
0017 }
0018
0019 CombinedMultiHitGenerator::~CombinedMultiHitGenerator() {}
0020
0021 void CombinedMultiHitGenerator::hitSets(const TrackingRegion& region,
0022 OrderedMultiHits& result,
0023 const edm::Event& ev,
0024 const edm::EventSetup& es) {
0025 edm::Handle<SeedingLayerSetsHits> hlayers;
0026 ev.getByToken(theSeedingLayerToken, hlayers);
0027 const SeedingLayerSetsHits& layers = *hlayers;
0028 if (layers.numberOfLayersInSet() != 3)
0029 throw cms::Exception("Configuration")
0030 << "CombinedMultiHitGenerator expects SeedingLayerSetsHits::numberOfLayersInSet() to be 3, got "
0031 << layers.numberOfLayersInSet();
0032
0033 theGenerator->initES(es);
0034
0035 std::vector<LayerTriplets::LayerSetAndLayers> trilayers = LayerTriplets::layers(layers);
0036 for (const auto& setAndLayers : trilayers) {
0037 theGenerator->hitSets(region, result, ev, es, setAndLayers.first, setAndLayers.second);
0038 }
0039 theLayerCache.clear();
0040 }