CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CosmicTrackSplitter.cc
Go to the documentation of this file.
5 
8 
14 
21 
29 
30 // added by me
31 
38 
46 
48 
51 
53 
56 
57 
61 
62 
63 #include <boost/regex.hpp>
64 
86 namespace reco { namespace modules {
88  public:
89  CosmicTrackSplitter(const edm::ParameterSet &iConfig) ;
90  virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) ;
91 
92  private:
96  size_t minimumHits_;
97 
103 
104  double dZcut_;
105  double dXYcut_;
106 
107  std::vector<uint32_t> detsToIgnore_;
108 
111 
112  TrackCandidate makeCandidate(const reco::Track &tk, std::vector<TrackingRecHit *>::iterator hitsBegin, std::vector<TrackingRecHit *>::iterator hitsEnd) ;
113 
114  }; // class
115 
116 
118  tracks_(iConfig.getParameter<edm::InputTag>("tracks")),
119  tjTag_(iConfig.getParameter<edm::InputTag>("tjTkAssociationMapTag") ),
120  minimumHits_(iConfig.getParameter<uint32_t>("minimumHits")),
121  replaceWithInactiveHits_(iConfig.getParameter<bool>("replaceWithInactiveHits")),
122  stripFrontInvalidHits_(iConfig.getParameter<bool>("stripFrontInvalidHits")),
123  stripBackInvalidHits_( iConfig.getParameter<bool>("stripBackInvalidHits") ),
124  stripAllInvalidHits_( iConfig.getParameter<bool>("stripAllInvalidHits") ),
125  excludePixelHits_( iConfig.getParameter<bool>("excludePixelHits") ),
126  dZcut_(iConfig.getParameter<double>("dzCut") ),
127  dXYcut_(iConfig.getParameter<double>("dxyCut") ),
128  detsToIgnore_( iConfig.getParameter<std::vector<uint32_t> >("detsToIgnore") )
129  {
130  // sanity check
132  throw cms::Exception("Configuration") << "Inconsistent Configuration: you can't set both 'stripAllInvalidHits' and 'replaceWithInactiveHits' to true\n";
133  }
134 
135  LogDebug("CosmicTrackSplitter") << "sanity check";
136 
137  // sort detids to ignore
138  std::sort(detsToIgnore_.begin(), detsToIgnore_.end());
139 
140  totalTracks_ = 0;
141 
142  // issue the produce<>
143  produces<TrackCandidateCollection>();
144  }
145 
146  void
148  {
149  LogDebug("CosmicTrackSplitter") << "IN THE SPLITTER!!!!!";
150 
151  // read with View, so we can read also a TrackRefVector
153  iEvent.getByLabel(tracks_, tracks);
154 
155  // also need trajectories ...
156  // Retrieve trajectories and tracks from the event
157  // -> merely skip if collection is empty
159  iEvent.getByLabel( tjTag_, m_TrajTracksMap );
160 
161  // read from EventSetup
164 
165  // prepare output collection
166  std::auto_ptr<TrackCandidateCollection> output(new TrackCandidateCollection());
167  output->reserve(tracks->size());
168 
169  // working area and tools
170  std::vector<TrackingRecHit *> hits;
171 
172  // Form pairs of trajectories and tracks
173  //ConstTrajTrackPairCollection trajTracks;
174  LogDebug("CosmicTrackSplitter") << "size of map: " << m_TrajTracksMap->size();
175  int HITTOSPLITFROM = 0;
176  for ( TrajTrackAssociationCollection::const_iterator iPair = m_TrajTracksMap->begin(); iPair != m_TrajTracksMap->end(); iPair++ ){
177  const Trajectory* trajFromMap = &(*(*iPair).key);
178  const reco::Track* trackFromMap = &(*(*iPair).val);
179 
180  // loop to find the hit to split from (by taking dot product of pT and transverse position
181  std::vector<TrajectoryMeasurement> measurements = trajFromMap->measurements();
182  int totalNumberOfHits = measurements.size();
183  int numberOfHits = 0;
184  double previousDotProduct = 0;
185  for (trackingRecHit_iterator ith = trackFromMap->recHitsBegin(), edh = trackFromMap->recHitsEnd(); ith != edh; ++ith) {
186 
187  GlobalVector stateMomentum = measurements[numberOfHits].forwardPredictedState().globalMomentum();
188  GlobalPoint statePosition = measurements[numberOfHits].forwardPredictedState().globalPosition();
189  double dotProduct = stateMomentum.x()*statePosition.x() + stateMomentum.y()*statePosition.y();
190  if ( dotProduct*previousDotProduct < 0 ){
191  //found hit to split from...
192  HITTOSPLITFROM = numberOfHits;
193  }
194 
195  previousDotProduct = dotProduct;
196  numberOfHits++;
197 
198  }
199  LogDebug("CosmicTrackSplitter") << "number of rechits: " << numberOfHits;
200 
201  // check if the trajectories and rechits are in reverse order...
202  trackingRecHit_iterator bIt = trackFromMap->recHitsBegin();
203  trackingRecHit_iterator fIt = trackFromMap->recHitsEnd() - 1;
204  const TrackingRecHit* bHit = bIt->get();
205  const TrackingRecHit* fHit = fIt->get();
206  // hit type valid = 0, missing = 1, inactive = 2, bad = 3
207  if( bHit->type() != 0 || bHit->isValid() != 1){
208  //loop over hits forwards until first Valid hit is found
210  for( ihit = trackFromMap->recHitsBegin();
211  ihit != trackFromMap->recHitsEnd(); ++ihit){
212  const TrackingRecHit* iHit = ihit->get();
213  if( iHit->type() == 0 && iHit->isValid() == 1){
214  bHit = iHit;
215  break;
216  }
217  }
218  }
219  DetId bdetid = bHit->geographicalId();
220  GlobalPoint bPosHit = theGeometry->idToDetUnit( bdetid)->surface().
221  toGlobal(bHit->localPosition());
222  if( fHit->type() != 0 || fHit->isValid() != 1){
223  //loop over hits backwards until first Valid hit is found
225  for( ihitf = trackFromMap->recHitsEnd()-1;
226  ihitf != trackFromMap->recHitsBegin(); --ihitf){
227  const TrackingRecHit* iHit = ihitf->get();
228  if( iHit->type() == 0 && iHit->isValid() == 1){
229  fHit = iHit;
230  break;
231  }
232  }
233  }
234  DetId fdetid = fHit->geographicalId();
235  GlobalPoint fPosHit = theGeometry->
236  idToDetUnit( fdetid )->surface().toGlobal(fHit->localPosition());
237  GlobalPoint bPosState = measurements[0].updatedState().globalPosition();
238  GlobalPoint fPosState = measurements[measurements.size()-1].
239  updatedState().globalPosition();
240  bool trajReversedFlag = false;
241  /*
242  DetId bdetid = bHit->geographicalId();
243  DetId fdetid = fHit->geographicalId();
244  GlobalPoint bPosHit = theGeometry->idToDetUnit( bdetid )->surface().toGlobal(bHit->localPosition());
245  GlobalPoint fPosHit = theGeometry->idToDetUnit( fdetid )->surface().toGlobal(fHit->localPosition());
246  GlobalPoint bPosState = measurements[0].updatedState().globalPosition();
247  GlobalPoint fPosState = measurements[measurements.size() - 1].updatedState().globalPosition();
248  bool trajReversedFlag = false;
249  */
250  if (( (bPosHit - bPosState).mag() > (bPosHit - fPosState).mag() ) && ( (fPosHit - fPosState).mag() > (fPosHit - bPosState).mag() ) ){
251  trajReversedFlag = true;
252  }
253  if (trajReversedFlag){ int temp = HITTOSPLITFROM; HITTOSPLITFROM = totalNumberOfHits - temp; }
254  }
255 
256  totalTracks_ = totalTracks_ + tracks->size();
257  // loop on tracks
258  for (std::vector<reco::Track>::const_iterator itt = tracks->begin(), edt = tracks->end(); itt != edt; ++itt) {
259  hits.clear(); // extra safety
260 
261  LogDebug("CosmicTrackSplitter") << "ntracks: " << tracks->size();
262 
263  // try to find distance of closest approach
264  math::XYZPoint refPoint = itt->referencePoint();
265  GlobalPoint v( itt->vx(), itt->vy(), itt->vz() );
266 
267  //checks on impact parameter
268  bool continueWithTrack = true;
269  if (fabs(v.z()) > dZcut_) continueWithTrack = false;
270  if (v.perp() > dXYcut_) continueWithTrack = false;
271  if (continueWithTrack == false) return;
272 
273  // LOOP TWICE, ONCE FOR TOP AND ONCE FOR BOTTOM
274  for (int i = 0; i < 2; ++i){
275  hits.clear(); // extra safety
276  LogDebug("CosmicTrackSplitter") << " loop on hits of track #" << (itt - tracks->begin());
277  int usedHitCtr = 0;
278  int hitCtr = 0;
279  for (trackingRecHit_iterator ith = itt->recHitsBegin(), edh = itt->recHitsEnd(); ith != edh; ++ith) {
280  //hitCtr++;
281  const TrackingRecHit * hit = ith->get(); // ith is an iterator on edm::Ref to rechit
282  LogDebug("CosmicTrackSplitter") << " hit number " << (ith - itt->recHitsBegin());
283  // let's look at valid hits
284  if (hit->isValid()) {
285  LogDebug("CosmicTrackSplitter") << " valid, detid = " << hit->geographicalId().rawId();
286  DetId detid = hit->geographicalId();
287 
288  if (detid.det() == DetId::Tracker) { // check for tracker hits
289  LogDebug("CosmicTrackSplitter") << " valid, tracker ";
290  bool verdict = false;
291 
292  //trying to get the global position of the hit
293  //const GeomDetUnit* geomDetUnit = theGeometry->idToDetUnit( detid ).;
294 
295  const GlobalPoint pos = theGeometry->idToDetUnit( detid )->surface().toGlobal(hit->localPosition());
296  LogDebug("CosmicTrackSplitter") << "hit pos: " << pos << ", dca pos: " << v;
297 
298  // top half
299  if ((i == 0)&&(hitCtr < HITTOSPLITFROM)){
300  verdict = true;
301  LogDebug("CosmicTrackSplitter") << "tophalf";
302  }
303  // bottom half
304  if ((i == 1)&&(hitCtr >= HITTOSPLITFROM)){
305  verdict = true;
306  LogDebug("CosmicTrackSplitter") << "bottomhalf";
307  }
308 
309  // if the hit is good, check again at module level
310  if ( verdict && std::binary_search(detsToIgnore_.begin(), detsToIgnore_.end(), detid.rawId())) {
311  verdict = false;
312  }
313 
314  // if hit is good check to make sure that we are keeping pixel hits
315  if ( excludePixelHits_){
316  if ((detid.det() == DetId::Tracker)&&((detid.subdetId() == 1)||(detid.subdetId() == 2))) { // check for pixel hits
317  verdict = false;
318  }
319  }
320 
321  LogDebug("CosmicTrackSplitter") << " verdict after module list: " << (verdict ? "ok" : "no");
322  if (verdict == true) {
323  // just copy the hit
324  hits.push_back(hit->clone());
325  usedHitCtr++;
326  }
327  else {
328  // still, if replaceWithInactiveHits is true we have to put a new hit
330  hits.push_back(new InvalidTrackingRecHit(detid, TrackingRecHit::inactive));
331  }
332  }
333  }
334  else { // just copy non tracker hits
335  hits.push_back(hit->clone());
336  }
337  }
338  else {
339  if (!stripAllInvalidHits_) {
340  hits.push_back(hit->clone());
341  }
342  } // is valid hit
343  LogDebug("CosmicTrackSplitter") << " end of hit " << (ith - itt->recHitsBegin());
344  hitCtr++;
345  } // loop on hits
346  LogDebug("CosmicTrackSplitter") << " end of loop on hits of track #" << (itt - tracks->begin());
347 
348  std::vector<TrackingRecHit *>::iterator begin = hits.begin(), end = hits.end();
349 
350  LogDebug("CosmicTrackSplitter") << " selected " << hits.size() << " hits ";
351 
352  // strip invalid hits at the beginning
354  while ( (begin != end) && ( (*begin)->isValid() == false ) ) ++begin;
355  }
356 
357  LogDebug("CosmicTrackSplitter") << " after front stripping we have " << (end - begin) << " hits ";
358 
359  // strip invalid hits at the end
360  if (stripBackInvalidHits_ && (begin != end)) {
361  --end;
362  while ( (begin != end) && ( (*end)->isValid() == false ) ) --end;
363  ++end;
364  }
365 
366  LogDebug("CosmicTrackSplitter") << " after back stripping we have " << (end - begin) << " hits ";
367 
368  // if we still have some hits
369  //if ((end - begin) >= int(minimumHits_)) {
370  if ( usedHitCtr >= int(minimumHits_)) {
371  output->push_back( makeCandidate ( *itt, begin, end ) );
372  LogDebug("CosmicTrackSplitter") << "we made a candidate of " << hits.size() << " hits!";
373  }
374  // now delete the hits not used by the candidate
375  for (begin = hits.begin(), end = hits.end(); begin != end; ++begin) {
376  if (*begin) delete *begin;
377  }
378  LogDebug("CosmicTrackSplitter") << "loop: " << i << " has " << usedHitCtr << " active hits and " << hits.size() << " total hits...";
379  hits.clear();
380  } // loop twice for top and bottom
381  } // loop on tracks
382  LogDebug("CosmicTrackSplitter") << "totalTracks_ = " << totalTracks_;
383  iEvent.put(output);
384  }
385 
387  CosmicTrackSplitter::makeCandidate(const reco::Track &tk, std::vector<TrackingRecHit *>::iterator hitsBegin, std::vector<TrackingRecHit *>::iterator hitsEnd) {
388 
389  LogDebug("CosmicTrackSplitter") << "Making a candidate!";
390 
391  TrajectoryStateTransform transform;
394  if ( pdir == anyDirection ) throw cms::Exception("UnimplementedFeature") << "Cannot work with tracks that have 'anyDirecton' \n";
395  //if ( (pdir == alongMomentum) == ( tk.p() >= tk.outerP() ) ) {
396  if ( (pdir == alongMomentum) == ( (tk.outerPosition()-tk.innerPosition()).Dot(tk.momentum()) >= 0 ) ) {
397  // use inner state
398  TrajectoryStateOnSurface originalTsosIn(transform.innerStateOnSurface(tk, *theGeometry, &*theMagField));
399  state = transform.persistentState( originalTsosIn, DetId(tk.innerDetId()) );
400  } else {
401  // use outer state
402  TrajectoryStateOnSurface originalTsosOut(transform.outerStateOnSurface(tk, *theGeometry, &*theMagField));
403  state = transform.persistentState( originalTsosOut, DetId(tk.outerDetId()) );
404  }
405 
406  TrajectorySeed seed(*state, TrackCandidate::RecHitContainer(), pdir);
407 
409  ownHits.reserve(hitsEnd - hitsBegin);
410  for ( ; hitsBegin != hitsEnd; ++hitsBegin) { ownHits.push_back( *hitsBegin ); }
411 
412  TrackCandidate cand(ownHits, seed, *state, tk.seedRef());
413  delete state;
414 
415  LogDebug("CosmicTrackSplitter") << " dumping the hits now: ";
416  for (TrackCandidate::range hitR = cand.recHits(); hitR.first != hitR.second; ++hitR.first) {
417  LogTrace("CosmicTrackSplitter") << " hit detid = " << hitR.first->geographicalId().rawId() <<
418  ", type = " << typeid(*hitR.first).name();
419  }
420 
421  return cand;
422  }
423 
424 }} //namespaces
425 
426 
427 // ========= MODULE DEF ==============
430 
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:150
TrackCandidate makeCandidate(const reco::Track &tk, std::vector< TrackingRecHit * >::iterator hitsBegin, std::vector< TrackingRecHit * >::iterator hitsEnd)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< TrackCandidate > TrackCandidateCollection
Type type() const
edm::ESHandle< TrackerGeometry > theGeometry
PropagationDirection
std::pair< const_iterator, const_iterator > range
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:47
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:42
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
DataContainer const & measurements() const
Definition: Trajectory.h:169
void push_back(D *&d)
Definition: OwnVector.h:290
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:59
edm::ESHandle< MagneticField > theMagField
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field) const
PTrajectoryStateOnDet * persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid) const
#define end
Definition: vmac.h:38
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
#define LogTrace(id)
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalVector
vector in glovbal coordinate system
Definition: Vector3D.h:28
virtual TrackingRecHit * clone() const =0
Definition: DetId.h:20
edm::RefToBase< TrajectorySeed > seedRef() const
Definition: Track.h:112
tuple tracks
Definition: testEve_cfg.py:39
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
const T & get() const
Definition: EventSetup.h:55
bool isValid() const
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field) const
char state
Definition: procUtils.cc:75
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:18
#define begin
Definition: vmac.h:31
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:105
CosmicTrackSplitter(const edm::ParameterSet &iConfig)
const reco::PFCandidateRefVector & tracks_
DetId geographicalId() const
Detector det() const
get the detector field from this detid
Definition: DetId.h:37
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
virtual LocalPoint localPosition() const =0
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:61
mathSSE::Vec4< T > v
void reserve(size_t)
Definition: OwnVector.h:284
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65