CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrajectorySeedProducer.h
Go to the documentation of this file.
1 #ifndef FastSimulation_Tracking_TrajectorySeedProducer_h
2 #define FastSimulation_Tracking_TrajectorySeedProducer_h
3 
11 
18 
23 
24 #include <memory>
25 #include <vector>
26 #include <sstream>
27 
28 
29 class MagneticField;
30 class MagneticFieldMap;
31 class TrackerGeometry;
34 
37 {
38  private:
40 
45  std::shared_ptr<PropagatorWithMaterial> thePropagator;
46 
50  std::unique_ptr<SeedCreator> seedCreator;
51  unsigned int minLayersCrossed;
52 
53  std::vector<std::vector<TrackingLayer>> seedingLayers;
54  double originRadius;
55  double ptMin;
57  double nSigmaZ;
58 
64  // tokens
67  public:
69 
71  {
72  }
73 
74  virtual void produce(edm::Event& e, const edm::EventSetup& es);
75 
77 
83 
91  inline bool passHitTuplesCuts(
92  const SeedingNode<TrackingLayer>& seedingNode,
93  const std::vector<TrajectorySeedHitCandidate>& trackerRecHits,
94  const std::vector<int>& hitIndicesInTree,
95  const TrajectorySeedHitCandidate& currentTrackerHit
96  ) const
97  {
98  switch (seedingNode.getDepth())
99  {
100  case 0:
101  {
102  return true;
103  /* example for 1 hits
104  const TrajectorySeedHitCandidate& hit1 = currentTrackerHit;
105  return pass1HitsCuts(hit1,trackingAlgorithmId);
106  */
107  }
108 
109  case 1:
110  {
111  const SeedingNode<TrackingLayer>* parentNode = &seedingNode;
112  parentNode = parentNode->getParent();
113  const TrajectorySeedHitCandidate& hit1 = trackerRecHits[hitIndicesInTree[parentNode->getIndex()]];
114  const TrajectorySeedHitCandidate& hit2 = currentTrackerHit;
115 
116  return pass2HitsCuts(hit1,hit2);
117  }
118  case 2:
119  {
120  return true;
121  /* example for 3 hits
122  const SeedingNode<LayerSpec>* parentNode = &seedingNode;
123  parentNode = parentNode->getParent();
124  const TrajectorySeedHitCandidate& hit2 = trackerRecHits[hitIndicesInTree[parentNode->getIndex()]];
125  parentNode = parentNode->getParent();
126  const TrajectorySeedHitCandidate& hit1 = trackerRecHits[hitIndicesInTree[parentNode->getIndex()]];
127  const TrajectorySeedHitCandidate& hit3 = currentTrackerHit;
128  return pass3HitsCuts(hit1,hit2,hit3,trackingAlgorithmId);
129  */
130  }
131  }
132  return true;
133  }
134 
135  bool pass2HitsCuts(const TrajectorySeedHitCandidate& hit1, const TrajectorySeedHitCandidate& hit2) const;
136 
138 
145  virtual std::vector<unsigned int> iterateHits(
146  unsigned int start,
147  const std::vector<TrajectorySeedHitCandidate>& trackerRecHits,
148  std::vector<int> hitIndicesInTree,
149  bool processSkippedHits
150  ) const;
151 
152  inline bool isHitOnLayer(const TrajectorySeedHitCandidate& trackerRecHit, const TrackingLayer& layer) const
153  {
154  return layer==trackerRecHit.getTrackingLayer();
155  }
156 
160  const GlobalPoint& gpos1,
161  const GlobalPoint& gpos2,
162  double error,
163  bool forward
164  ) const;
165 
169  const GlobalPoint& gpos1,
170  const GlobalPoint& gpos2,
171  double error,
172  bool forward
173  ) const;
174 
176 
183  bool testWithRegions(const TrajectorySeedHitCandidate & innerHit,const TrajectorySeedHitCandidate & outerHit) const;
185  const std::vector<TrajectorySeedHitCandidate>& trackerRecHits,
186  std::vector<int>& hitIndicesInTree,
187  const SeedingNode<TrackingLayer>* node,
188  unsigned int trackerHit
189  ) const;
190 
191  // typedef std::vector<TrackingRegion* > Regions;
192  typedef std::vector<std::unique_ptr<TrackingRegion> > Regions;
194  //TrackingRegionProducer* theRegionProducer;
195 
196  std::unique_ptr<TrackingRegionProducer> theRegionProducer;
200 
201 
202 };
203 
204 #endif
bool passHitTuplesCuts(const SeedingNode< TrackingLayer > &seedingNode, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, const std::vector< int > &hitIndicesInTree, const TrajectorySeedHitCandidate &currentTrackerHit) const
method checks if a SimTrack fulfills the quality requirements.
edm::EDGetTokenT< FastTrackerRecHitCombinationCollection > recHitCombinationsToken
tuple start
Check for commandline option errors.
Definition: dqm_diff.py:58
const MagneticField * magneticField
bool compatibleWithBeamSpot(const GlobalPoint &gpos1, const GlobalPoint &gpos2, double error, bool forward) const
const SeedingNode< TrackingLayer > * insertHit(const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > &hitIndicesInTree, const SeedingNode< TrackingLayer > *node, unsigned int trackerHit) const
TrajectorySeedProducer(const edm::ParameterSet &conf)
virtual void produce(edm::Event &e, const edm::EventSetup &es)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::unique_ptr< SeedCreator > seedCreator
edm::EDGetTokenT< std::vector< bool > > hitMasksToken
unsigned int getDepth() const
Definition: SeedingTree.h:84
const reco::BeamSpot * beamSpot
std::vector< std::vector< TrackingLayer > > seedingLayers
const SeedingNode * getParent() const
Definition: SeedingTree.h:112
std::unique_ptr< TrackingRegionProducer > theRegionProducer
const TrackerTopology * trackerTopology
const edm::EventSetup * es_
const reco::VertexCollection * primaryVertices
bool isHitOnLayer(const TrajectorySeedHitCandidate &trackerRecHit, const TrackingLayer &layer) const
bool pass2HitsCuts(const TrajectorySeedHitCandidate &hit1, const TrajectorySeedHitCandidate &hit2) const
SeedingTree< TrackingLayer > _seedingTree
virtual std::vector< unsigned int > iterateHits(unsigned int start, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > hitIndicesInTree, bool processSkippedHits) const
method tries to insert all hits into the tree structure.
const MeasurementTrackerEvent * measurementTrackerEvent
std::vector< std::unique_ptr< TrackingRegion > > Regions
edm::EDGetTokenT< MeasurementTrackerEvent > measurementTrackerEventToken
const TrackerGeometry * trackerGeometry
std::shared_ptr< PropagatorWithMaterial > thePropagator
unsigned int getIndex() const
Definition: SeedingTree.h:107
const TrackingLayer & getTrackingLayer() const
bool compatibleWithPrimaryVertex(const GlobalPoint &gpos1, const GlobalPoint &gpos2, double error, bool forward) const
const MagneticFieldMap * magneticFieldMap
bool testWithRegions(const TrajectorySeedHitCandidate &innerHit, const TrajectorySeedHitCandidate &outerHit) const
method inserts hit into the tree structure at an empty position.