CMS 3D CMS Logo

SeedFinderSelector.cc
Go to the documentation of this file.
2 
3 // framework
6 
7 // track reco
21 // data formats
23 
25  : trackingRegion_(nullptr)
26  , eventSetup_(nullptr)
27  , measurementTracker_(nullptr)
28  , measurementTrackerLabel_(cfg.getParameter<std::string>("measurementTracker"))
29 {
30  if(cfg.exists("pixelTripletGeneratorFactory"))
31  {
32  const edm::ParameterSet & tripletConfig = cfg.getParameter<edm::ParameterSet>("pixelTripletGeneratorFactory");
33  pixelTripletGenerator_.reset(HitTripletGeneratorFromPairAndLayersFactory::get()->create(tripletConfig.getParameter<std::string>("ComponentName"),tripletConfig,consumesCollector));
34  }
35 
36  if(cfg.exists("MultiHitGeneratorFactory"))
37  {
38  const edm::ParameterSet & tripletConfig = cfg.getParameter<edm::ParameterSet>("MultiHitGeneratorFactory");
39  multiHitGenerator_.reset(MultiHitGeneratorFromPairAndLayersFactory::get()->create(tripletConfig.getParameter<std::string>("ComponentName"),tripletConfig));
40  }
41 
42  if(cfg.exists("CAHitTripletGeneratorFactory"))
43  {
44  const edm::ParameterSet & tripletConfig = cfg.getParameter<edm::ParameterSet>("CAHitTripletGeneratorFactory");
45  CAHitTriplGenerator_ = std::make_unique<CAHitTripletGenerator>(tripletConfig,consumesCollector);
46  seedingLayers_ = std::make_unique<SeedingLayerSetsBuilder>(cfg, consumesCollector, edm::InputTag("fastTrackerRecHits")); //calling the new FastSim specific constructor to make SeedingLayerSetsHits pointer for triplet iterations
47  layerPairs_ = cfg.getParameter<std::vector<unsigned>>("layerPairs"); //allowed layer pairs for CA triplets
48  }
49 
50  if(cfg.exists("CAHitQuadrupletGeneratorFactory"))
51  {
52  const edm::ParameterSet & quadrupletConfig = cfg.getParameter<edm::ParameterSet>("CAHitQuadrupletGeneratorFactory");
53  CAHitQuadGenerator_ = std::make_unique<CAHitQuadrupletGenerator>(quadrupletConfig, consumesCollector);
54  seedingLayers_ = std::make_unique<SeedingLayerSetsBuilder>(cfg, consumesCollector, edm::InputTag("fastTrackerRecHits")); //calling the new FastSim specific constructor to make SeedingLayerSetsHits pointer for quadruplet iterations
55  layerPairs_ = cfg.getParameter<std::vector<unsigned>>("layerPairs"); //allowed layer pairs for CA quadruplets
56  }
57 
59  {
60  throw cms::Exception("FastSimTracking") << "It is forbidden to specify together 'pixelTripletGeneratorFactory', 'CAHitTripletGeneratorFactory' and 'MultiHitGeneratorFactory' in configuration of SeedFinderSelection";
61  }
63  {
64  throw cms::Exception("FastSimTracking") << "It is forbidden to specify 'CAHitQuadrupletGeneratorFactory' together with 'pixelTripletGeneratorFactory', 'CAHitTripletGeneratorFactory' or 'MultiHitGeneratorFactory' in configuration of SeedFinderSelection";
65  }
66 }
67 
68 
70 
72 {
73  eventSetup_ = &es;
74 
75  edm::ESHandle<MeasurementTracker> measurementTrackerHandle;
76  es.get<CkfComponentsRecord>().get(measurementTrackerLabel_, measurementTrackerHandle);
78  measurementTracker_ = &(*measurementTrackerHandle);
79 
81  {
82  multiHitGenerator_->initES(es);
83  }
84 
85  //for CA triplet iterations
87  seedingLayer = seedingLayers_->makeSeedingLayerSetsHitsforFastSim(ev, es);
88  CAHitTriplGenerator_->initEvent(ev,es);
89  }
90  //for CA quadruplet iterations
92  seedingLayer = seedingLayers_->makeSeedingLayerSetsHitsforFastSim(ev, es);
93  CAHitQuadGenerator_->initEvent(ev,es);
94  }
95 }
96 
97 
98 bool SeedFinderSelector::pass(const std::vector<const FastTrackerRecHit *>& hits) const
99 {
101  {
102  throw cms::Exception("FastSimTracking") << "ERROR: event not initialized";
103  }
104  if(!trackingRegion_)
105  {
106  throw cms::Exception("FastSimTracking") << "ERROR: trackingRegion not set";
107  }
108 
109 
110  // check the inner 2 hits
111  if(hits.size() < 2)
112  {
113  throw cms::Exception("FastSimTracking") << "SeedFinderSelector::pass requires at least 2 hits";
114  }
115  const DetLayer * firstLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[0]->det()->geographicalId());
116  const DetLayer * secondLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[1]->det()->geographicalId());
117 
118  std::vector<BaseTrackerRecHit const *> firstHits{hits[0]};
119  std::vector<BaseTrackerRecHit const *> secondHits{hits[1]};
120 
121  const RecHitsSortedInPhi fhm(firstHits, trackingRegion_->origin(), firstLayer);
122  const RecHitsSortedInPhi shm(secondHits, trackingRegion_->origin(), secondLayer);
123 
124  HitDoublets result(fhm,shm);
125  HitPairGeneratorFromLayerPair::doublets(*trackingRegion_,*firstLayer,*secondLayer,fhm,shm,*eventSetup_,0,result);
126 
127  if(result.empty())
128  {
129  return false;
130  }
131 
132  // check the inner 3 hits
134  {
135  if(hits.size() < 3)
136  {
137  throw cms::Exception("FastSimTracking") << "For the given configuration, SeedFinderSelector::pass requires at least 3 hits";
138  }
139  const DetLayer * thirdLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[2]->det()->geographicalId());
140  std::vector<const DetLayer *> thirdLayerDetLayer(1,thirdLayer);
141  std::vector<BaseTrackerRecHit const *> thirdHits{hits[2]};
142  const RecHitsSortedInPhi thm(thirdHits,trackingRegion_->origin(), thirdLayer);
143  const RecHitsSortedInPhi * thmp =&thm;
144 
146  {
147  OrderedHitTriplets tripletresult;
148  pixelTripletGenerator_->hitTriplets(*trackingRegion_,tripletresult,*eventSetup_,result,&thmp,thirdLayerDetLayer,1);
149  return !tripletresult.empty();
150  }
151  else if(multiHitGenerator_)
152  {
153  OrderedMultiHits tripletresult;
154  multiHitGenerator_->hitTriplets(*trackingRegion_,tripletresult,*eventSetup_,result,&thmp,thirdLayerDetLayer,1);
155  return !tripletresult.empty();
156  }
157  //new for Phase1
158  else if(CAHitTriplGenerator_)
159  {
160  if(!seedingLayer)
161  throw cms::Exception("FastSimTracking") << "ERROR: SeedingLayers pointer not set for CATripletGenerator";
162 
164  //constructing IntermediateHitDoublets to be passed onto CAHitTripletGenerator::hitNtuplets()
165  IntermediateHitDoublets ihd(&layers);
166  const TrackingRegion& tr_ = *trackingRegion_;
167  auto filler = ihd.beginRegion(&tr_);
168 
169  //Forming doublets for CA triplets from the allowed layer pair combinations:(0,1),(1,2)
170  std::string hitPair[3];
171  hitPair[0] = Layer_name(hits[0]);
172  hitPair[1] = Layer_name(hits[1]);
173  hitPair[2] = Layer_name(hits[2]);
174 
175  //extracting the DetLayer of the hits
176  const DetLayer * fLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[0]->det()->geographicalId());
177  const DetLayer * sLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[1]->det()->geographicalId());
178  const DetLayer * tLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[2]->det()->geographicalId());
179 
180  //converting FastTrackerRecHit hits to BaseTrackerRecHit
181  std::vector<BaseTrackerRecHit const *> fHits{hits[0]};
182  std::vector<BaseTrackerRecHit const *> sHits{hits[1]};
183  std::vector<BaseTrackerRecHit const *> tHits{hits[2]};
184 
185  //forming the SeedingLayerSet for the hit doublets
186  SeedingLayerSetsHits::SeedingLayerSet pairCandidate1, pairCandidate2;
189  for(const auto p : layerPairs_){
190  pairCandidate = ls.slice(p,p+2);
191  if(p==0 && hitPair[0] == pairCandidate[0].name() && hitPair[1] == pairCandidate[1].name())
192  pairCandidate1=pairCandidate;
193  if(p==1 && hitPair[1] == pairCandidate[0].name() && hitPair[2] == pairCandidate[1].name())
194  pairCandidate2=pairCandidate;
195  }
196  }
197 
198  //Important: hits of the layer to be added to LayerHitMapCache
199  auto& layerCache = filler.layerHitMapCache();
200 
201  //doublets for CA triplets from the allowed layer pair combinations:(0,1),(1,2) and storing in filler
202  const RecHitsSortedInPhi & firsthm = *layerCache.add(pairCandidate1[0], std::make_unique<RecHitsSortedInPhi>(fHits, trackingRegion_->origin(), fLayer));
203  const RecHitsSortedInPhi & secondhm = *layerCache.add(pairCandidate1[1], std::make_unique<RecHitsSortedInPhi>(sHits, trackingRegion_->origin(), sLayer));
204  HitDoublets res1(firsthm,secondhm);
205  HitPairGeneratorFromLayerPair::doublets(*trackingRegion_,*fLayer,*sLayer,firsthm,secondhm,*eventSetup_,0,res1);
206  filler.addDoublets(pairCandidate1, std::move(res1));
207  const RecHitsSortedInPhi & thirdhm = *layerCache.add(pairCandidate2[1], std::make_unique<RecHitsSortedInPhi>(tHits, trackingRegion_->origin(), tLayer));
208  HitDoublets res2(secondhm,thirdhm);
209  HitPairGeneratorFromLayerPair::doublets(*trackingRegion_,*sLayer,*tLayer,secondhm,thirdhm,*eventSetup_,0,res2);
210  filler.addDoublets(pairCandidate2, std::move(res2));
211 
212  std::vector<OrderedHitSeeds> tripletresult;
213  tripletresult.resize(ihd.regionSize());
214  for(auto& ntuplet : tripletresult)
215  ntuplet.reserve(3);
216  CAHitTriplGenerator_->hitNtuplets(ihd,tripletresult,*eventSetup_,*seedingLayer); //calling the function from the class, modifies tripletresult
217  return !tripletresult.empty();
218  }
219  }
220  //new for Phase1
222  {
223  if(hits.size() < 4)
224  {
225  throw cms::Exception("FastSimTracking") << "For the given configuration, SeedFinderSelector::pass requires at least 4 hits";
226  }
227 
228  if(!seedingLayer)
229  throw cms::Exception("FastSimTracking") << "ERROR: SeedingLayers pointer not set for CAHitQuadrupletGenerator";
230 
232  //constructing IntermediateHitDoublets to be passed onto CAHitQuadrupletGenerator::hitNtuplets()
233  IntermediateHitDoublets ihd(&layers);
234  const TrackingRegion& tr_ = *trackingRegion_;
235  auto filler = ihd.beginRegion(&tr_);
236 
237  //Forming doublets for CA quadruplets from the allowed layer pair combinations:(0,1),(1,2),(2,3)
238  std::string hitPair[4];
239  hitPair[0] = Layer_name(hits[0]);
240  hitPair[1] = Layer_name(hits[1]);
241  hitPair[2] = Layer_name(hits[2]);
242  hitPair[3] = Layer_name(hits[3]);
243 
244  //extracting the DetLayer of the hits
245  const DetLayer * fLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[0]->det()->geographicalId());
246  const DetLayer * sLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[1]->det()->geographicalId());
247  const DetLayer * tLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[2]->det()->geographicalId());
248  const DetLayer * frLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[3]->det()->geographicalId());
249 
250  //converting FastTrackerRecHit hits to BaseTrackerRecHit
251  std::vector<BaseTrackerRecHit const *> fHits{hits[0]};
252  std::vector<BaseTrackerRecHit const *> sHits{hits[1]};
253  std::vector<BaseTrackerRecHit const *> tHits{hits[2]};
254  std::vector<BaseTrackerRecHit const *> frHits{hits[3]};
255 
256  //forming the SeedingLayerSet for the hit doublets
257  SeedingLayerSetsHits::SeedingLayerSet pairCandidate1, pairCandidate2, pairCandidate3;
260  for(const auto p : layerPairs_){
261  pairCandidate = ls.slice(p,p+2);
262  if(p==0 && hitPair[0] == pairCandidate[0].name() && hitPair[1] == pairCandidate[1].name())
263  pairCandidate1=pairCandidate;
264  if(p==1 && hitPair[1] == pairCandidate[0].name() && hitPair[2] == pairCandidate[1].name())
265  pairCandidate2=pairCandidate;
266  if(p==2 && hitPair[2] == pairCandidate[0].name() && hitPair[3] == pairCandidate[1].name())
267  pairCandidate3=pairCandidate;
268  }
269  }
270 
271  //Important: hits of the layer to be added to LayerHitMapCache
272  auto& layerCache = filler.layerHitMapCache();
273 
274  //doublets for CA quadruplets from the allowed layer pair combinations:(0,1),(1,2),(2,3) and storing in filler
275  const RecHitsSortedInPhi & firsthm = *layerCache.add(pairCandidate1[0], std::make_unique<RecHitsSortedInPhi>(fHits, trackingRegion_->origin(), fLayer));
276  const RecHitsSortedInPhi & secondhm = *layerCache.add(pairCandidate1[1], std::make_unique<RecHitsSortedInPhi>(sHits, trackingRegion_->origin(), sLayer));
277  HitDoublets res1(firsthm,secondhm);
278  HitPairGeneratorFromLayerPair::doublets(*trackingRegion_,*fLayer,*sLayer,firsthm,secondhm,*eventSetup_,0,res1);
279  filler.addDoublets(pairCandidate1, std::move(res1));
280  const RecHitsSortedInPhi & thirdhm = *layerCache.add(pairCandidate2[1], std::make_unique<RecHitsSortedInPhi>(tHits, trackingRegion_->origin(), tLayer));
281  HitDoublets res2(secondhm,thirdhm);
282  HitPairGeneratorFromLayerPair::doublets(*trackingRegion_,*sLayer,*tLayer,secondhm,thirdhm,*eventSetup_,0,res2);
283  filler.addDoublets(pairCandidate2, std::move(res2));
284  const RecHitsSortedInPhi & fourthhm = *layerCache.add(pairCandidate3[1], std::make_unique<RecHitsSortedInPhi>(frHits, trackingRegion_->origin(), frLayer));
285  HitDoublets res3(thirdhm, fourthhm);
286  HitPairGeneratorFromLayerPair::doublets(*trackingRegion_,*tLayer,*frLayer,thirdhm,fourthhm,*eventSetup_,0,res3);
287  filler.addDoublets(pairCandidate3, std::move(res3));
288 
289  std::vector<OrderedHitSeeds> quadrupletresult;
290  quadrupletresult.resize(ihd.regionSize());
291  for(auto& ntuplet : quadrupletresult)
292  ntuplet.reserve(4);
293  CAHitQuadGenerator_->hitNtuplets(ihd,quadrupletresult,*eventSetup_,*seedingLayer); //calling the function from the class, modifies quadrupletresult
294  return !quadrupletresult.empty();
295  }
296 
297  return true;
298 
299 }
300 
301 //new for Phase1
303 {
304  std::string name="UNKWN";
305  const TrackerTopology* const tTopo = trackerTopology.product();
306  int idLayer = 0;
307  if( (hit->det()->geographicalId()).subdetId() == PixelSubdetector::PixelBarrel){
308  idLayer = tTopo->pxbLayer(hit->det()->geographicalId());
309  name = "BPix"+std::to_string(idLayer);
310  }
311  else if ((hit->det()->geographicalId()).subdetId() == PixelSubdetector::PixelEndcap){
312  idLayer = tTopo->pxfDisk(hit->det()->geographicalId());
313  if(tTopo->pxfSide(hit->det()->geographicalId())==1)
314  name = "FPix"+std::to_string(idLayer)+"_neg";
315  else
316  name = "FPix"+std::to_string(idLayer)+"_pos";
317  }
318  return name;
319 }
T getParameter(std::string const &) const
bool pass(const std::vector< const FastTrackerRecHit * > &hits) const
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
std::unique_ptr< CAHitQuadrupletGenerator > CAHitQuadGenerator_
const edm::EventSetup * eventSetup_
GlobalPoint const & origin() const
def create(alignables, pedeDump, additionalData, outputFile, config)
std::unique_ptr< SeedingLayerSetsHits > seedingLayer
unsigned int pxfDisk(const DetId &id) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool ev
SeedFinderSelector(const edm::ParameterSet &cfg, edm::ConsumesCollector &&consumesCollector)
RegionFiller beginRegion(const TrackingRegion *region)
#define nullptr
std::unique_ptr< HitTripletGeneratorFromPairAndLayers > pixelTripletGenerator_
std::string Layer_name(const FastTrackerRecHit *hit) const
const GeomDet * det() const
std::vector< unsigned > layerPairs_
const TrackingRegion * trackingRegion_
const DetLayer * detLayer(const DetId &id) const
obsolete method. Use idToLayer() instead.
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
HitDoublets doublets(const TrackingRegion &reg, const edm::Event &ev, const edm::EventSetup &es, Layers layers)
edm::ESHandle< TrackerTopology > trackerTopology
unsigned int pxbLayer(const DetId &id) const
def ls(path, rec=False)
Definition: eostools.py:348
const T & get() const
Definition: EventSetup.h:55
void initEvent(const edm::Event &ev, const edm::EventSetup &es)
std::unique_ptr< SeedingLayerSetsBuilder > seedingLayers_
const MeasurementTracker * measurementTracker_
unsigned int pxfSide(const DetId &id) const
const GeometricSearchTracker * geometricSearchTracker() const
const std::string measurementTrackerLabel_
T const * product() const
Definition: ESHandle.h:86
std::unique_ptr< CAHitTripletGenerator > CAHitTriplGenerator_
def move(src, dest)
Definition: eostools.py:510
T get(const Candidate &c)
Definition: component.h:55
std::unique_ptr< MultiHitGeneratorFromPairAndLayers > multiHitGenerator_