CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackCandidateProducer.cc
Go to the documentation of this file.
1 // system
2 #include <memory>
3 #include <vector>
4 #include <map>
5 
6 // framework
14 
15 // data format
24 
25 // geometry / magnetic field / propagation
35 
36 // fastsim
41 
43 {
44 public:
45 
46  explicit TrackCandidateProducer(const edm::ParameterSet& conf);
47 
48  virtual void produce(edm::Event& e, const edm::EventSetup& es) override;
49 
50 private:
51 
52  // tokens & labels
58 
59  // other data
61  bool splitHits;
64 };
65 
67  : hitSplitter()
68 {
69  // produces
70  produces<TrackCandidateCollection>();
71 
72  // consumes
73  if (conf.exists("hitMasks")){
74  hitMasksToken = consumes<std::vector<bool> >(conf.getParameter<edm::InputTag>("hitMasks"));
75  }
76  seedToken = consumes<edm::View<TrajectorySeed> >(conf.getParameter<edm::InputTag>("src"));
77  recHitCombinationsToken = consumes<FastTrackerRecHitCombinationCollection>(conf.getParameter<edm::InputTag>("recHitCombinations"));
78  simTrackToken = consumes<edm::SimTrackContainer>(conf.getParameter<edm::InputTag>("simTracks"));
79 
80  // other parameters
81  maxSeedMatchEstimator = conf.getUntrackedParameter<double>("maxSeedMatchEstimator",0);
82  rejectOverlaps = conf.getParameter<bool>("OverlapCleaning");
83  splitHits = conf.getParameter<bool>("SplitHits");
84  propagatorLabel = conf.getParameter<std::string>("propagator");
85 }
86 
87 void
89 
90  // get records
92  es.get<IdealMagneticFieldRecord>().get(magneticField);
93 
94  edm::ESHandle<TrackerGeometry> trackerGeometry;
95  es.get<TrackerDigiGeometryRecord>().get(trackerGeometry);
96 
97  edm::ESHandle<TrackerTopology> trackerTopology;
98  es.get<TrackerTopologyRcd>().get(trackerTopology);
99 
101  es.get<TrackingComponentsRecord>().get(propagatorLabel,propagator);
102 
103  // get products
105  e.getByToken(seedToken,seeds);
106 
108  e.getByToken(recHitCombinationsToken, recHitCombinations);
109 
111  e.getByToken(simTrackToken,simTracks);
112 
113  const std::vector<bool> * hitMasks = 0;
115  edm::Handle<std::vector<bool> > hitMasksHandle;
116  e.getByToken(hitMasksToken,hitMasksHandle);
117  hitMasks = &(*hitMasksHandle);
118  }
119 
120  // output collection
121  std::unique_ptr<TrackCandidateCollection> output(new TrackCandidateCollection);
122 
123  // loop over the seeds
124  for (unsigned seedIndex = 0; seedIndex < seeds->size(); ++seedIndex){
125 
126  const TrajectorySeed seed = (*seeds)[seedIndex];
127  std::vector<int32_t> recHitCombinationIndices;
128 
129  // match hitless seeds to simTracks and find corresponding recHitCombination
130  if(seed.nHits()==0){
131  recHitCombinationIndices = SeedMatcher::matchRecHitCombinations(
132  seed,
133  *recHitCombinations,
134  *simTracks,
136  *propagator,
137  *magneticField,
138  *trackerGeometry);
139  }
140  // for normal seeds, retrieve the corresponding recHitCombination from the seed hits
141  else
142  {
144  recHitCombinationIndices.push_back(icomb);
145  }
146 
147  // loop over the matched recHitCombinations
148  for(auto icomb : recHitCombinationIndices)
149  {
150  if(icomb < 0 || unsigned(icomb) >= recHitCombinations->size())
151  {
152  throw cms::Exception("TrackCandidateProducer") << " found seed with recHitCombination out or range: " << icomb << std::endl;
153  }
154  const FastTrackerRecHitCombination & recHitCombination = (*recHitCombinations)[icomb];
155 
156  // container for select hits
157  std::vector<const FastTrackerRecHit *> selectedRecHits;
158 
159  // add the seed hits
160  TrajectorySeed::range seedHitRange = seed.recHits();//Hits in a seed
161  for (TrajectorySeed::const_iterator ihit = seedHitRange.first; ihit != seedHitRange.second; ++ihit)
162  {
163  selectedRecHits.push_back(static_cast<const FastTrackerRecHit*>(&*ihit));
164  }
165 
166  // prepare to skip seed hits
167  const FastTrackerRecHit * lastHitToSkip = 0;
168  if(selectedRecHits.size() > 0)
169  {
170  lastHitToSkip = selectedRecHits.back();
171  }
172 
173  // inOut or outIn tracking ?
174  bool hitsAlongMomentum = (seed.direction()== alongMomentum);
175 
176  // add hits from combination to hit selection
177  for (unsigned hitIndex = hitsAlongMomentum ? 0 : recHitCombination.size() - 1;
178  hitIndex < recHitCombination.size();
179  hitsAlongMomentum ? ++hitIndex : --hitIndex)
180  {
181 
182  const FastTrackerRecHit * selectedRecHit = recHitCombination[hitIndex].get();
183 
184  // skip seed hits
185  if(lastHitToSkip)
186  {
187  if(lastHitToSkip->sameId(selectedRecHit))
188  {
189  lastHitToSkip=0;
190  }
191  continue;
192  }
193 
194  // apply hit masking
195  if(hitMasks && fastTrackingUtilities::hitIsMasked(selectedRecHit,hitMasks))
196  {
197  continue;
198  }
199 
200 
201  // if overlap rejection is not switched on, accept all hits
202  // always accept the first hit
203  // also accept a hit if it is not on the layer of the previous hit
204  if( ! rejectOverlaps
205  || selectedRecHits.size() == 0
206  || ( TrackingLayer::createFromDetId(selectedRecHits.back()->geographicalId(),*trackerTopology.product())
207  != TrackingLayer::createFromDetId(selectedRecHit->geographicalId(),*trackerTopology.product())))
208  {
209  selectedRecHits.push_back(selectedRecHit);
210  }
211  // else:
212  // overlap rejection is switched on
213  // the hit is on the same layer as the previous hit
214  // accept the one with smallest error
215  else if ( fastTrackingUtilities::hitLocalError(selectedRecHit)
216  < fastTrackingUtilities::hitLocalError(selectedRecHits.back()) )
217  {
218  selectedRecHits.back() = selectedRecHit;
219  }
220  }
221 
222  // split hits / store copies for the track candidate
223  edm::OwnVector<TrackingRecHit> hitsForTrackCandidate;
224  for ( unsigned index = 0; index<selectedRecHits.size(); ++index )
225  {
226  if(splitHits)
227  {
228  // add split hits to splitSelectedRecHits
229  hitSplitter.split(*selectedRecHits[index],hitsForTrackCandidate,hitsAlongMomentum);
230  }
231  else
232  {
233  hitsForTrackCandidate.push_back(selectedRecHits[index]->clone());
234  }
235  }
236 
237  // set the recHitCombinationIndex
238  fastTrackingUtilities::setRecHitCombinationIndex(hitsForTrackCandidate,icomb);
239 
240  // create track candidate state
241  // 1. get seed state (defined on the surface of the most outer hit)
242  DetId seedDetId(seed.startingState().detId());
243  const GeomDet* gdet = trackerGeometry->idToDet(seedDetId);
244  TrajectoryStateOnSurface seedTSOS = trajectoryStateTransform::transientState(seed.startingState(), &(gdet->surface()),magneticField.product());
245  // 2. backPropagate the seedState to the surfuce of the most inner hit
246  const GeomDet* initialLayer = trackerGeometry->idToDet(hitsForTrackCandidate.front().geographicalId());
247  const TrajectoryStateOnSurface initialTSOS = propagator->propagate(seedTSOS,initialLayer->surface()) ;
248  // 3. check validity and transform
249  if (!initialTSOS.isValid()) continue;
250  PTrajectoryStateOnDet PTSOD = trajectoryStateTransform::persistentState(initialTSOS,hitsForTrackCandidate.front().geographicalId().rawId());
251 
252  // add track candidate to output collection
253  output->push_back(TrackCandidate(hitsForTrackCandidate,seed,PTSOD,edm::RefToBase<TrajectorySeed>(seeds,seedIndex)));
254  }
255  }
256 
257  // Save the track candidates
258  e.put(std::move(output));
259 
260 }
261 
PropagationDirection direction() const
void setRecHitCombinationIndex(edm::OwnVector< T > &recHits, int32_t icomb)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
static std::vector< int > matchRecHitCombinations(const TrajectorySeed &seed, const FastTrackerRecHitCombinationCollection &recHitCombinationCollection, const std::vector< SimTrack > &simTrackCollection, double maxMatchEstimator, const Propagator &propagator, const MagneticField &magneticField, const TrackerGeometry &trackerGeometry)
Definition: SeedMatcher.cc:14
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
bool sameId(const FastTrackerRecHit *other, size_t i=0, size_t j=0) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
TrackCandidateProducer(const edm::ParameterSet &conf)
static TrackingLayer createFromDetId(const DetId &detId, const TrackerTopology &trackerTopology)
std::vector< TrackCandidate > TrackCandidateCollection
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool hitIsMasked(const FastTrackerRecHit *hit, const std::vector< bool > *hitMasks)
edm::EDGetTokenT< edm::View< TrajectorySeed > > seedToken
FastTrackerRecHitSplitter hitSplitter
edm::EDGetTokenT< FastTrackerRecHitCombinationCollection > recHitCombinationsToken
edm::EDGetTokenT< std::vector< bool > > hitMasksToken
virtual void produce(edm::Event &e, const edm::EventSetup &es) override
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
void push_back(D *&d)
Definition: OwnVector.h:290
const SurfaceType & surface() const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
recHitContainer::const_iterator const_iterator
def move
Definition: eostools.py:510
std::pair< const_iterator, const_iterator > range
unsigned int detId() const
int32_t getRecHitCombinationIndex(const T &object)
Definition: DetId.h:18
std::vector< FastTrackerRecHitRef > FastTrackerRecHitCombination
PTrajectoryStateOnDet const & startingState() const
void split(const FastTrackerRecHit &hitIn, edm::OwnVector< TrackingRecHit > &hitsOut, bool alongMomentum) const
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
range recHits() const
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
unsigned int nHits() const
double hitLocalError(const TrackingRecHit *hit)
bool isUninitialized() const
Definition: EDGetToken.h:73
edm::EDGetTokenT< edm::SimTrackContainer > simTrackToken
DetId geographicalId() const
reference front()
Definition: OwnVector.h:429