CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CosmicTrackFinder.cc
Go to the documentation of this file.
1 // Package: RecoTracker/SingleTrackPattern
2 // Class: CosmicTrackFinder
3 // Original Author: Michele Pioppi-INFN perugia
4 #include <memory>
5 #include <string>
6 
23 
25 
26 
27 namespace cms
28 {
29 
31  cosmicTrajectoryBuilder_(conf) ,
32  crackTrajectoryBuilder_(conf)
33  {
34  geometry=conf.getUntrackedParameter<std::string>("GeometricStructure","STANDARD");
35  useHitsSplitting_=conf.getParameter<bool>("useHitsSplitting");
36  matchedrecHitsToken_ = consumes<SiStripMatchedRecHit2DCollection>(
37  conf.getParameter<edm::InputTag>("matchedRecHits"));
38  rphirecHitsToken_ = consumes<SiStripRecHit2DCollection>(
39  conf.getParameter<edm::InputTag>("rphirecHits"));
40  stereorecHitsToken_ = consumes<SiStripRecHit2DCollection>(
41  conf.getParameter<edm::InputTag>("stereorecHits"));
42  pixelRecHitsToken_ = consumes<SiPixelRecHitCollection>(
43  conf.getParameter<edm::InputTag>("pixelRecHits"));
44  seedToken_ = consumes<TrajectorySeedCollection>(
45  conf.getParameter<edm::InputTag>("cosmicSeeds"));
46 
47  produces<TrackCandidateCollection>();
48  }
49 
50 
51  // Virtual destructor needed.
53 
54  // Functions that gets called by framework every event
56  {
57  using namespace std;
58 
59  // retrieve seeds
62 
63  //retrieve PixelRecHits
64  static const SiPixelRecHitCollection s_empty;
67  if (geometry!="MTCC" && (geometry!="CRACK" )) {
68  if( e.getByToken(pixelRecHitsToken_, pixelHits)) {
69  pixelHitCollection = pixelHits.product();
70  } else {
71  Labels l;
73  edm::LogWarning("CosmicTrackFinder")
74  << "Collection SiPixelRecHitCollection with InputTag "
75  << l.module
76  << " cannot be found, using empty collection of same type.";
77  }
78  }
79 
80 
81 
82 
83  //retrieve StripRecHits
85  e.getByToken(matchedrecHitsToken_, matchedrecHits);
87  e.getByToken(rphirecHitsToken_, rphirecHits);
89  e.getByToken(stereorecHitsToken_, stereorecHits);
90 
91  // Step B: create empty output collection
92  std::auto_ptr<TrackCandidateCollection> output(new TrackCandidateCollection);
93 
96  edm::LogVerbatim("CosmicTrackFinder") << "========== Cosmic Track Finder Info ==========";
97  edm::LogVerbatim("CosmicTrackFinder") << " Numbers of Seeds " << (*seed).size();
98  if((*seed).size()>0){
99 
100  std::vector<Trajectory> trajoutput;
101 
102 
103  const TransientTrackingRecHitBuilder * hitBuilder = nullptr;
104  if(geometry!="CRACK" ) {
106  *stereorecHits,
107  *rphirecHits,
108  *matchedrecHits,
110  es,
111  e,
112  trajoutput);
113  hitBuilder = cosmicTrajectoryBuilder_.hitBuilder();
114  } else {
116  *stereorecHits,
117  *rphirecHits,
118  *matchedrecHits,
120  es,
121  e,
122  trajoutput);
123  hitBuilder = crackTrajectoryBuilder_.hitBuilder();
124  }
125  assert(hitBuilder);
126  Traj2TrackHits t2t(hitBuilder,true);
127 
128  edm::LogVerbatim("CosmicTrackFinder") << " Numbers of Temp Trajectories " << trajoutput.size();
129  edm::LogVerbatim("CosmicTrackFinder") << "========== END Info ==========";
130  if(trajoutput.size()>0){
131 
132  // crazyness...
133  std::vector<Trajectory*> tmpTraj;
134  std::vector<Trajectory>::iterator itr;
135  for (itr=trajoutput.begin();itr!=trajoutput.end();itr++)tmpTraj.push_back(&(*itr));
136 
137  //The best track is selected
138  //FOR MTCC the criteria are:
139  //1)# of layers,2) # of Hits,3)Chi2
140  if (geometry=="MTCC") stable_sort(tmpTraj.begin(),tmpTraj.end(),CompareTrajLay());
141  else stable_sort(tmpTraj.begin(),tmpTraj.end(),CompareTrajChi());
142 
144 
145  const Trajectory theTraj = *(*tmpTraj.begin());
146  bool seedplus=(theTraj.seed().direction()==alongMomentum);
147 
148  // std::cout << "CosmicTrackFinder " <<"Reconstruction " << (seedplus ? "along" : "opposite to") << " momentum" << std::endl;
149  LogDebug("CosmicTrackFinder")<<"Reconstruction " << (seedplus ? "along" : "opposite to") << " momentum";
150 
151  /*
152  // === the convention is to save always final tracks with hits sorted *along* momentum
153  // --- this is NOT what was necessaraly happening before and not consistent with what done in standard CKF (this is a candidate not a track)
154  */
156  if(theTraj.direction() == alongMomentum) std::cout << "cosmic: along momentum... " << std::endl;
157  t2t(theTraj,recHits,useHitsSplitting_);
158  recHits.reverse(); // according to original code
159 
160  /*
161  Trajectory::RecHitContainer thits;
162  //it->recHitsV(thits);
163  theTraj.recHitsV(thits,useHitsSplitting_);
164  edm::OwnVector<TrackingRecHit> recHits;
165  recHits.reserve(thits.size());
166  // reverse hit order
167  for (Trajectory::RecHitContainer::const_iterator hitIt = thits.end()-1;
168  hitIt >= thits.begin(); hitIt--) {
169  recHits.push_back( (**hitIt).hit()->clone());
170  }
171  */
172 
173  TSOS firstState;
174  unsigned int firstId;
175 
176  // assume not along momentum....
177  firstState=theTraj.lastMeasurement().updatedState();
178  firstId = theTraj.lastMeasurement().recHitR().rawId();
179  //firstId = recHits.front().rawId();
180 
181  /*
182  cout << "firstState y, z: " << firstState.globalPosition().y()
183  << " , " << firstState.globalPosition().z() << endl;
184 
185  */
186 
188  TSOS startingState( firstState.localParameters(), LocalTrajectoryError(C),
189  firstState.surface(),
190  firstState.magneticField() );
191 
192  // protection againt invalid initial states
193  if (! firstState.isValid()) {
194  edm::LogWarning("CosmicTrackFinder") << "invalid innerState, will not make TrackCandidate";
195  edm::OrphanHandle<TrackCandidateCollection> rTrackCand = e.put(output);
196  return;
197  }
198 
199  // FIXME in case of slitting this can happen see CkfTrackCandidateMakerBase
200  if(firstId != recHits.front().rawId()){
201  edm::LogWarning("CosmicTrackFinder") <<"Mismatch in DetID of first hit: firstID= " <<firstId
202  << " DetId= " << recHits.front().geographicalId().rawId();
203  edm::OrphanHandle<TrackCandidateCollection> rTrackCand = e.put(output);
204  return;
205  }
206 
207  PTrajectoryStateOnDet const & state = trajectoryStateTransform::persistentState( startingState, firstId);
208 
209 
210  output->push_back(TrackCandidate(recHits,theTraj.seed(),state,theTraj.seedRef() ) );
211 
212  }
213 
214  }
215  edm::OrphanHandle<TrackCandidateCollection> rTrackCand = e.put(output);
216  }
217 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const LocalTrajectoryParameters & localParameters() const
edm::EDGetTokenT< SiStripRecHit2DCollection > stereorecHitsToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
ROOT::Math::SMatrixIdentity AlgebraicMatrixID
assert(m_qm.get())
std::vector< TrackCandidate > TrackCandidateCollection
CosmicTrajectoryBuilder cosmicTrajectoryBuilder_
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const MagneticField * magneticField() const
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const TransientTrackingRecHitBuilder * hitBuilder() const
void reverse()
Definition: OwnVector.h:159
const SurfaceType & surface() const
void run(const TrajectorySeedCollection &collseed, const SiStripRecHit2DCollection &collstereo, const SiStripRecHit2DCollection &collrphi, const SiStripMatchedRecHit2DCollection &collmatched, const SiPixelRecHitCollection &collpixel, const edm::EventSetup &es, edm::Event &e, std::vector< Trajectory > &trajoutput)
Runs the algorithm.
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
CRackTrajectoryBuilder crackTrajectoryBuilder_
edm::EDGetTokenT< SiPixelRecHitCollection > pixelRecHitsToken_
CosmicTrackFinder(const edm::ParameterSet &conf)
char const * module
Definition: ProductLabels.h:5
const TransientTrackingRecHitBuilder * hitBuilder() const
edm::EDGetTokenT< TrajectorySeedCollection > seedToken_
edm::EDGetTokenT< SiStripMatchedRecHit2DCollection > matchedrecHitsToken_
virtual void produce(edm::Event &e, const edm::EventSetup &c)
void run(const TrajectorySeedCollection &collseed, const SiStripRecHit2DCollection &collstereo, const SiStripRecHit2DCollection &collrphi, const SiStripMatchedRecHit2DCollection &collmatched, const SiPixelRecHitCollection &collpixel, const edm::EventSetup &es, edm::Event &e, std::vector< Trajectory > &trajoutput)
Runs the algorithm.
const T & get() const
Definition: EventSetup.h:56
edm::EDGetTokenT< SiStripRecHit2DCollection > rphirecHitsToken_
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
tuple cout
Definition: gather_cfg.py:145
DetId geographicalId() const
reference front()
Definition: OwnVector.h:429
id_type rawId() const