Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:11:23

0001 #include "FastSimulation/Tracking/interface/SeedFinderSelector.h"
0002 
0003 // framework
0004 #include "FWCore/Framework/interface/ConsumesCollector.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 
0008 // track reco
0009 #include "RecoTracker/TkHitPairs/interface/RecHitsSortedInPhi.h"
0010 #include "RecoTracker/TkHitPairs/interface/HitPairGeneratorFromLayerPair.h"
0011 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0012 #include "RecoTracker/TkSeedGenerator/interface/MultiHitGeneratorFromPairAndLayers.h"
0013 #include "RecoTracker/TkSeedGenerator/interface/MultiHitGeneratorFromPairAndLayersFactory.h"
0014 #include "RecoTracker/PixelSeeding/interface/HitTripletGeneratorFromPairAndLayers.h"
0015 #include "RecoTracker/PixelSeeding/interface/HitTripletGeneratorFromPairAndLayersFactory.h"
0016 #include "RecoTracker/PixelSeeding/interface/CAHitTripletGenerator.h"
0017 #include "RecoTracker/PixelSeeding/interface/CAHitQuadrupletGenerator.h"
0018 #include "RecoTracker/PixelSeeding/interface/OrderedHitSeeds.h"
0019 #include "RecoTracker/TkHitPairs/interface/IntermediateHitDoublets.h"
0020 #include "RecoTracker/TkMSParametrization/interface/MultipleScatteringParametrisationMaker.h"
0021 #include "RecoTracker/Record/interface/TrackerMultipleScatteringRecord.h"
0022 // data formats
0023 #include "DataFormats/TrackerRecHit2D/interface/FastTrackerRecHit.h"
0024 
0025 SeedFinderSelector::SeedFinderSelector(const edm::ParameterSet &cfg, edm::ConsumesCollector &&consumesCollector)
0026     : trackingRegion_(nullptr),
0027       eventSetup_(nullptr),
0028       measurementTracker_(nullptr),
0029       measurementTrackerLabel_(cfg.getParameter<std::string>("measurementTracker")),
0030       measurementTrackerESToken_(consumesCollector.esConsumes(edm::ESInputTag("", measurementTrackerLabel_))),
0031       trackerTopologyESToken_(consumesCollector.esConsumes()),
0032       fieldESToken_(consumesCollector.esConsumes()),
0033       msMakerESToken_(consumesCollector.esConsumes()) {
0034   if (cfg.exists("pixelTripletGeneratorFactory")) {
0035     const edm::ParameterSet &tripletConfig = cfg.getParameter<edm::ParameterSet>("pixelTripletGeneratorFactory");
0036     pixelTripletGenerator_ = HitTripletGeneratorFromPairAndLayersFactory::get()->create(
0037         tripletConfig.getParameter<std::string>("ComponentName"), tripletConfig, consumesCollector);
0038   }
0039 
0040   if (cfg.exists("MultiHitGeneratorFactory")) {
0041     const edm::ParameterSet &tripletConfig = cfg.getParameter<edm::ParameterSet>("MultiHitGeneratorFactory");
0042     multiHitGenerator_ = MultiHitGeneratorFromPairAndLayersFactory::get()->create(
0043         tripletConfig.getParameter<std::string>("ComponentName"), tripletConfig, consumesCollector);
0044   }
0045 
0046   if (cfg.exists("CAHitTripletGeneratorFactory")) {
0047     const edm::ParameterSet &tripletConfig = cfg.getParameter<edm::ParameterSet>("CAHitTripletGeneratorFactory");
0048     CAHitTriplGenerator_ = std::make_unique<CAHitTripletGenerator>(tripletConfig, consumesCollector);
0049     seedingLayers_ = std::make_unique<SeedingLayerSetsBuilder>(
0050         cfg,
0051         consumesCollector,
0052         //calling the new FastSim specific constructor to make SeedingLayerSetsHits pointer for triplet iterations
0053         edm::InputTag("fastTrackerRecHits"));
0054     layerPairs_ = cfg.getParameter<std::vector<unsigned>>("layerPairs");  //allowed layer pairs for CA triplets
0055   }
0056 
0057   if (cfg.exists("CAHitQuadrupletGeneratorFactory")) {
0058     const edm::ParameterSet &quadrupletConfig = cfg.getParameter<edm::ParameterSet>("CAHitQuadrupletGeneratorFactory");
0059     CAHitQuadGenerator_ = std::make_unique<CAHitQuadrupletGenerator>(quadrupletConfig, consumesCollector);
0060     //calling the new FastSim specific constructor to make SeedingLayerSetsHits pointer for quadruplet iterations
0061     seedingLayers_ =
0062         std::make_unique<SeedingLayerSetsBuilder>(cfg, consumesCollector, edm::InputTag("fastTrackerRecHits"));
0063     layerPairs_ = cfg.getParameter<std::vector<unsigned>>("layerPairs");  //allowed layer pairs for CA quadruplets
0064   }
0065 
0066   if ((pixelTripletGenerator_ && multiHitGenerator_) || (CAHitTriplGenerator_ && pixelTripletGenerator_) ||
0067       (CAHitTriplGenerator_ && multiHitGenerator_)) {
0068     throw cms::Exception("FastSimTracking")
0069         << "It is forbidden to specify together 'pixelTripletGeneratorFactory', 'CAHitTripletGeneratorFactory' and "
0070            "'MultiHitGeneratorFactory' in configuration of SeedFinderSelection";
0071   }
0072   if ((pixelTripletGenerator_ && CAHitQuadGenerator_) || (CAHitTriplGenerator_ && CAHitQuadGenerator_) ||
0073       (multiHitGenerator_ && CAHitQuadGenerator_)) {
0074     throw cms::Exception("FastSimTracking")
0075         << "It is forbidden to specify 'CAHitQuadrupletGeneratorFactory' together with 'pixelTripletGeneratorFactory', "
0076            "'CAHitTripletGeneratorFactory' or 'MultiHitGeneratorFactory' in configuration of SeedFinderSelection";
0077   }
0078 }
0079 
0080 SeedFinderSelector::~SeedFinderSelector() { ; }
0081 
0082 void SeedFinderSelector::initEvent(const edm::Event &ev, const edm::EventSetup &es) {
0083   eventSetup_ = &es;
0084 
0085   measurementTracker_ = &es.getData(measurementTrackerESToken_);
0086   trackerTopology_ = &es.getData(trackerTopologyESToken_);
0087   field_ = &es.getData(fieldESToken_);
0088   msmaker_ = &es.getData(msMakerESToken_);
0089 
0090   if (multiHitGenerator_) {
0091     multiHitGenerator_->initES(es);
0092   }
0093 
0094   //for CA triplet iterations
0095   if (CAHitTriplGenerator_) {
0096     seedingLayer = seedingLayers_->makeSeedingLayerSetsHitsforFastSim(ev, es);
0097     seedingLayerIds = seedingLayers_->layers();
0098     CAHitTriplGenerator_->initEvent(ev, es);
0099   }
0100   //for CA quadruplet iterations
0101   if (CAHitQuadGenerator_) {
0102     seedingLayer = seedingLayers_->makeSeedingLayerSetsHitsforFastSim(ev, es);
0103     seedingLayerIds = seedingLayers_->layers();
0104     CAHitQuadGenerator_->initEvent(ev, es);
0105   }
0106 }
0107 
0108 bool SeedFinderSelector::pass(const std::vector<const FastTrackerRecHit *> &hits) const {
0109   if (!measurementTracker_ || !eventSetup_) {
0110     throw cms::Exception("FastSimTracking") << "ERROR: event not initialized";
0111   }
0112   if (!trackingRegion_) {
0113     throw cms::Exception("FastSimTracking") << "ERROR: trackingRegion not set";
0114   }
0115 
0116   // check the inner 2 hits
0117   if (hits.size() < 2) {
0118     throw cms::Exception("FastSimTracking") << "SeedFinderSelector::pass requires at least 2 hits";
0119   }
0120   const DetLayer *firstLayer =
0121       measurementTracker_->geometricSearchTracker()->detLayer(hits[0]->det()->geographicalId());
0122   const DetLayer *secondLayer =
0123       measurementTracker_->geometricSearchTracker()->detLayer(hits[1]->det()->geographicalId());
0124 
0125   std::vector<BaseTrackerRecHit const *> firstHits{hits[0]};
0126   std::vector<BaseTrackerRecHit const *> secondHits{hits[1]};
0127 
0128   const RecHitsSortedInPhi fhm(firstHits, trackingRegion_->origin(), firstLayer);
0129   const RecHitsSortedInPhi shm(secondHits, trackingRegion_->origin(), secondLayer);
0130 
0131   HitDoublets result(fhm, shm);
0132   HitPairGeneratorFromLayerPair::doublets(
0133       *trackingRegion_, *firstLayer, *secondLayer, fhm, shm, *field_, *msmaker_, 0, result);
0134 
0135   if (result.empty()) {
0136     return false;
0137   }
0138 
0139   // check the inner 3 hits
0140   if (pixelTripletGenerator_ || multiHitGenerator_ || CAHitTriplGenerator_) {
0141     if (hits.size() < 3) {
0142       throw cms::Exception("FastSimTracking")
0143           << "For the given configuration, SeedFinderSelector::pass requires at least 3 hits";
0144     }
0145     const DetLayer *thirdLayer =
0146         measurementTracker_->geometricSearchTracker()->detLayer(hits[2]->det()->geographicalId());
0147     std::vector<const DetLayer *> thirdLayerDetLayer(1, thirdLayer);
0148     std::vector<BaseTrackerRecHit const *> thirdHits{hits[2]};
0149     const RecHitsSortedInPhi thm(thirdHits, trackingRegion_->origin(), thirdLayer);
0150     const RecHitsSortedInPhi *thmp = &thm;
0151 
0152     if (pixelTripletGenerator_) {
0153       OrderedHitTriplets tripletresult;
0154       pixelTripletGenerator_->hitTriplets(
0155           *trackingRegion_, tripletresult, *eventSetup_, result, &thmp, thirdLayerDetLayer, 1);
0156       return !tripletresult.empty();
0157     } else if (multiHitGenerator_) {
0158       OrderedMultiHits tripletresult;
0159       multiHitGenerator_->hitTriplets(*trackingRegion_, tripletresult, result, &thmp, thirdLayerDetLayer, 1);
0160       return !tripletresult.empty();
0161     }
0162     //new for Phase1
0163     else if (CAHitTriplGenerator_) {
0164       if (!seedingLayer)
0165         throw cms::Exception("FastSimTracking") << "ERROR: SeedingLayers pointer not set for CATripletGenerator";
0166 
0167       SeedingLayerSetsHits &layers = *seedingLayer;
0168       //constructing IntermediateHitDoublets to be passed onto CAHitTripletGenerator::hitNtuplets()
0169       IntermediateHitDoublets ihd(&layers);
0170       const TrackingRegion &tr_ = *trackingRegion_;
0171       auto filler = ihd.beginRegion(&tr_);
0172 
0173       //forming the SeedingLayerId of the hits
0174       std::array<SeedingLayerSetsBuilder::SeedingLayerId, 3> hitPair;
0175       hitPair[0] = Layer_tuple(hits[0]);
0176       hitPair[1] = Layer_tuple(hits[1]);
0177       hitPair[2] = Layer_tuple(hits[2]);
0178 
0179       //extracting the DetLayer of the hits
0180       const DetLayer *fLayer =
0181           measurementTracker_->geometricSearchTracker()->detLayer(hits[0]->det()->geographicalId());
0182       const DetLayer *sLayer =
0183           measurementTracker_->geometricSearchTracker()->detLayer(hits[1]->det()->geographicalId());
0184       const DetLayer *tLayer =
0185           measurementTracker_->geometricSearchTracker()->detLayer(hits[2]->det()->geographicalId());
0186 
0187       //converting FastTrackerRecHit hits to BaseTrackerRecHit
0188       std::vector<BaseTrackerRecHit const *> fHits{hits[0]};
0189       std::vector<BaseTrackerRecHit const *> sHits{hits[1]};
0190       std::vector<BaseTrackerRecHit const *> tHits{hits[2]};
0191 
0192       //forming the SeedingLayerSet for the hit doublets
0193       SeedingLayerSetsHits::SeedingLayerSet pairCandidate1, pairCandidate2;
0194       for (SeedingLayerSetsHits::SeedingLayerSet ls : *seedingLayer) {
0195         SeedingLayerSetsHits::SeedingLayerSet pairCandidate;
0196         for (const auto p : layerPairs_) {
0197           pairCandidate = ls.slice(p, p + 2);
0198           if (p == 0 && hitPair[0] == seedingLayerIds[pairCandidate[0].index()] &&
0199               hitPair[1] == seedingLayerIds[pairCandidate[1].index()])
0200             pairCandidate1 = pairCandidate;
0201           if (p == 1 && hitPair[1] == seedingLayerIds[pairCandidate[0].index()] &&
0202               hitPair[2] == seedingLayerIds[pairCandidate[1].index()])
0203             pairCandidate2 = pairCandidate;
0204         }
0205       }
0206 
0207       //Important: hits of the layer to be added to LayerHitMapCache
0208       auto &layerCache = filler.layerHitMapCache();
0209 
0210       //doublets for CA triplets from the allowed layer pair combinations:(0,1),(1,2) and storing in filler
0211       const RecHitsSortedInPhi &firsthm = *layerCache.add(
0212           pairCandidate1[0], std::make_unique<RecHitsSortedInPhi>(fHits, trackingRegion_->origin(), fLayer));
0213       const RecHitsSortedInPhi &secondhm = *layerCache.add(
0214           pairCandidate1[1], std::make_unique<RecHitsSortedInPhi>(sHits, trackingRegion_->origin(), sLayer));
0215       HitDoublets res1(firsthm, secondhm);
0216       HitPairGeneratorFromLayerPair::doublets(
0217           *trackingRegion_, *fLayer, *sLayer, firsthm, secondhm, *field_, *msmaker_, 0, res1);
0218       filler.addDoublets(pairCandidate1, std::move(res1));
0219       const RecHitsSortedInPhi &thirdhm = *layerCache.add(
0220           pairCandidate2[1], std::make_unique<RecHitsSortedInPhi>(tHits, trackingRegion_->origin(), tLayer));
0221       HitDoublets res2(secondhm, thirdhm);
0222       HitPairGeneratorFromLayerPair::doublets(
0223           *trackingRegion_, *sLayer, *tLayer, secondhm, thirdhm, *field_, *msmaker_, 0, res2);
0224       filler.addDoublets(pairCandidate2, std::move(res2));
0225 
0226       std::vector<OrderedHitSeeds> tripletresult;
0227       tripletresult.resize(ihd.regionSize());
0228       for (auto &ntuplet : tripletresult)
0229         ntuplet.reserve(3);
0230       //calling the function from the class, modifies tripletresult
0231       CAHitTriplGenerator_->hitNtuplets(ihd, tripletresult, *seedingLayer);
0232       return !tripletresult[0].empty();
0233     }
0234   }
0235   //new for Phase1
0236   if (CAHitQuadGenerator_) {
0237     if (hits.size() < 4) {
0238       throw cms::Exception("FastSimTracking")
0239           << "For the given configuration, SeedFinderSelector::pass requires at least 4 hits";
0240     }
0241 
0242     if (!seedingLayer)
0243       throw cms::Exception("FastSimTracking") << "ERROR: SeedingLayers pointer not set for CAHitQuadrupletGenerator";
0244 
0245     SeedingLayerSetsHits &layers = *seedingLayer;
0246     //constructing IntermediateHitDoublets to be passed onto CAHitQuadrupletGenerator::hitNtuplets()
0247     IntermediateHitDoublets ihd(&layers);
0248     const TrackingRegion &tr_ = *trackingRegion_;
0249     auto filler = ihd.beginRegion(&tr_);
0250 
0251     //forming the SeedingLayerId of the hits
0252     std::array<SeedingLayerSetsBuilder::SeedingLayerId, 4> hitPair;
0253     hitPair[0] = Layer_tuple(hits[0]);
0254     hitPair[1] = Layer_tuple(hits[1]);
0255     hitPair[2] = Layer_tuple(hits[2]);
0256     hitPair[3] = Layer_tuple(hits[3]);
0257 
0258     //extracting the DetLayer of the hits
0259     const DetLayer *fLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[0]->det()->geographicalId());
0260     const DetLayer *sLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[1]->det()->geographicalId());
0261     const DetLayer *tLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[2]->det()->geographicalId());
0262     const DetLayer *frLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[3]->det()->geographicalId());
0263 
0264     //converting FastTrackerRecHit hits to BaseTrackerRecHit
0265     std::vector<BaseTrackerRecHit const *> fHits{hits[0]};
0266     std::vector<BaseTrackerRecHit const *> sHits{hits[1]};
0267     std::vector<BaseTrackerRecHit const *> tHits{hits[2]};
0268     std::vector<BaseTrackerRecHit const *> frHits{hits[3]};
0269 
0270     //forming the SeedingLayerSet for the hit doublets
0271     SeedingLayerSetsHits::SeedingLayerSet pairCandidate1, pairCandidate2, pairCandidate3;
0272     for (SeedingLayerSetsHits::SeedingLayerSet ls : *seedingLayer) {
0273       SeedingLayerSetsHits::SeedingLayerSet pairCandidate;
0274       for (const auto p : layerPairs_) {
0275         pairCandidate = ls.slice(p, p + 2);
0276         if (p == 0 && hitPair[0] == seedingLayerIds[pairCandidate[0].index()] &&
0277             hitPair[1] == seedingLayerIds[pairCandidate[1].index()])
0278           pairCandidate1 = pairCandidate;
0279         if (p == 1 && hitPair[1] == seedingLayerIds[pairCandidate[0].index()] &&
0280             hitPair[2] == seedingLayerIds[pairCandidate[1].index()])
0281           pairCandidate2 = pairCandidate;
0282         if (p == 2 && hitPair[2] == seedingLayerIds[pairCandidate[0].index()] &&
0283             hitPair[3] == seedingLayerIds[pairCandidate[1].index()])
0284           pairCandidate3 = pairCandidate;
0285       }
0286     }
0287 
0288     //Important: hits of the layer to be added to LayerHitMapCache
0289     auto &layerCache = filler.layerHitMapCache();
0290 
0291     //doublets for CA quadruplets from the allowed layer pair combinations:(0,1),(1,2),(2,3) and storing in filler
0292     const RecHitsSortedInPhi &firsthm = *layerCache.add(
0293         pairCandidate1[0], std::make_unique<RecHitsSortedInPhi>(fHits, trackingRegion_->origin(), fLayer));
0294     const RecHitsSortedInPhi &secondhm = *layerCache.add(
0295         pairCandidate1[1], std::make_unique<RecHitsSortedInPhi>(sHits, trackingRegion_->origin(), sLayer));
0296     HitDoublets res1(firsthm, secondhm);
0297     HitPairGeneratorFromLayerPair::doublets(
0298         *trackingRegion_, *fLayer, *sLayer, firsthm, secondhm, *field_, *msmaker_, 0, res1);
0299     filler.addDoublets(pairCandidate1, std::move(res1));
0300     const RecHitsSortedInPhi &thirdhm = *layerCache.add(
0301         pairCandidate2[1], std::make_unique<RecHitsSortedInPhi>(tHits, trackingRegion_->origin(), tLayer));
0302     HitDoublets res2(secondhm, thirdhm);
0303     HitPairGeneratorFromLayerPair::doublets(
0304         *trackingRegion_, *sLayer, *tLayer, secondhm, thirdhm, *field_, *msmaker_, 0, res2);
0305     filler.addDoublets(pairCandidate2, std::move(res2));
0306     const RecHitsSortedInPhi &fourthhm = *layerCache.add(
0307         pairCandidate3[1], std::make_unique<RecHitsSortedInPhi>(frHits, trackingRegion_->origin(), frLayer));
0308     HitDoublets res3(thirdhm, fourthhm);
0309     HitPairGeneratorFromLayerPair::doublets(
0310         *trackingRegion_, *tLayer, *frLayer, thirdhm, fourthhm, *field_, *msmaker_, 0, res3);
0311     filler.addDoublets(pairCandidate3, std::move(res3));
0312 
0313     std::vector<OrderedHitSeeds> quadrupletresult;
0314     quadrupletresult.resize(ihd.regionSize());
0315     for (auto &ntuplet : quadrupletresult)
0316       ntuplet.reserve(4);
0317     //calling the function from the class, modifies quadrupletresult
0318     CAHitQuadGenerator_->hitNtuplets(ihd, quadrupletresult, *seedingLayer);
0319     return !quadrupletresult[0].empty();
0320   }
0321 
0322   return true;
0323 }
0324 
0325 //new for Phase1
0326 SeedingLayerSetsBuilder::SeedingLayerId SeedFinderSelector::Layer_tuple(const FastTrackerRecHit *hit) const {
0327   GeomDetEnumerators::SubDetector subdet = GeomDetEnumerators::invalidDet;
0328   TrackerDetSide side = TrackerDetSide::Barrel;
0329   int idLayer = 0;
0330 
0331   if ((hit->det()->geographicalId()).subdetId() == PixelSubdetector::PixelBarrel) {
0332     subdet = GeomDetEnumerators::PixelBarrel;
0333     side = TrackerDetSide::Barrel;
0334     idLayer = trackerTopology_->pxbLayer(hit->det()->geographicalId());
0335   } else if ((hit->det()->geographicalId()).subdetId() == PixelSubdetector::PixelEndcap) {
0336     subdet = GeomDetEnumerators::PixelEndcap;
0337     idLayer = trackerTopology_->pxfDisk(hit->det()->geographicalId());
0338     if (trackerTopology_->pxfSide(hit->det()->geographicalId()) == 1) {
0339       side = TrackerDetSide::NegEndcap;
0340     } else {
0341       side = TrackerDetSide::PosEndcap;
0342     }
0343   }
0344   return std::make_tuple(subdet, side, idLayer);
0345 }