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 
10 
12 
16 
20 
23  //#include "DataFormats/BeamSpot/interface/BeamSpot.h"
24  //#include "DataFormats/TrackReco/interface/TrackFwd.h"
25 
26 #include <memory>
27 #include <vector>
28 #include <sstream>
29 
30 
31 class MagneticField;
32 class MagneticFieldMap;
33 class TrackerGeometry;
35 
38 {
39  private:
41 
46 
47  std::shared_ptr<PropagatorWithMaterial> thePropagator;
48 
52 
53  unsigned int minLayersCrossed;
54 
55  std::vector<std::vector<TrackingLayer>> seedingLayers;
56 
57  double originRadius;
58  double ptMin;
60  double nSigmaZ;
61 
66  // tokens
72  std::vector<edm::EDGetTokenT<std::vector<unsigned int> > > skipSimTrackIdTokens;
73 
74  public:
75 
77 
79  {
80  }
81 
82  virtual void beginRun(edm::Run const& run, const edm::EventSetup & es);
83  virtual void produce(edm::Event& e, const edm::EventSetup& es);
84 
86 
91  virtual bool passSimTrackQualityCuts(const SimTrack& theSimTrack, const SimVertex& theSimVertex) const;
92 
94 
101  inline bool passHitTuplesCuts(
102  const SeedingNode<TrackingLayer>& seedingNode,
103  const std::vector<TrajectorySeedHitCandidate>& trackerRecHits,
104  const std::vector<int>& hitIndicesInTree,
105  const TrajectorySeedHitCandidate& currentTrackerHit
106  ) const
107  {
108  switch (seedingNode.getDepth())
109  {
110  case 0:
111  {
112  return true;
113  /* example for 1 hits
114  const TrajectorySeedHitCandidate& hit1 = currentTrackerHit;
115  return pass1HitsCuts(hit1,trackingAlgorithmId);
116  */
117  }
118 
119  case 1:
120  {
121  const SeedingNode<TrackingLayer>* parentNode = &seedingNode;
122  parentNode = parentNode->getParent();
123  const TrajectorySeedHitCandidate& hit1 = trackerRecHits[hitIndicesInTree[parentNode->getIndex()]];
124  const TrajectorySeedHitCandidate& hit2 = currentTrackerHit;
125 
126  return pass2HitsCuts(hit1,hit2);
127  }
128  case 2:
129  {
130  return true;
131  /* example for 3 hits
132  const SeedingNode<LayerSpec>* parentNode = &seedingNode;
133  parentNode = parentNode->getParent();
134  const TrajectorySeedHitCandidate& hit2 = trackerRecHits[hitIndicesInTree[parentNode->getIndex()]];
135  parentNode = parentNode->getParent();
136  const TrajectorySeedHitCandidate& hit1 = trackerRecHits[hitIndicesInTree[parentNode->getIndex()]];
137  const TrajectorySeedHitCandidate& hit3 = currentTrackerHit;
138  return pass3HitsCuts(hit1,hit2,hit3,trackingAlgorithmId);
139  */
140  }
141  }
142  return true;
143  }
144 
145  bool pass2HitsCuts(const TrajectorySeedHitCandidate& hit1, const TrajectorySeedHitCandidate& hit2) const;
146 
148 
155  virtual std::vector<unsigned int> iterateHits(
156  unsigned int start,
157  const std::vector<TrajectorySeedHitCandidate>& trackerRecHits,
158  std::vector<int> hitIndicesInTree,
159  bool processSkippedHits
160  ) const;
161 
162  inline bool isHitOnLayer(const TrajectorySeedHitCandidate& trackerRecHit, const TrackingLayer& layer) const
163  {
164  return layer==trackerRecHit.getTrackingLayer();
165  }
166 
170  const GlobalPoint& gpos1,
171  const GlobalPoint& gpos2,
172  double error,
173  bool forward
174  ) const;
175 
179  const GlobalPoint& gpos1,
180  const GlobalPoint& gpos2,
181  double error,
182  bool forward
183  ) const;
184 
186 
194  const std::vector<TrajectorySeedHitCandidate>& trackerRecHits,
195  std::vector<int>& hitIndicesInTree,
197  unsigned int trackerHit
198  ) const;
199 
200 
201 };
202 
203 #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 TrajectorySeedHitCandidate fulfills the quality requirements.
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
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
edm::EDGetTokenT< SiTrackerGSMatchedRecHit2DCollection > recHitToken
tuple node
Definition: Node.py:50
unsigned int getDepth() const
Definition: SeedingTree.h:84
const reco::BeamSpot * beamSpot
std::vector< edm::EDGetTokenT< std::vector< unsigned int > > > skipSimTrackIdTokens
const SeedingNode * getParent() const
Definition: SeedingTree.h:112
const TrackerTopology * trackerTopology
const reco::VertexCollection * primaryVertices
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.
virtual void beginRun(edm::Run const &run, const edm::EventSetup &es)
edm::EDGetTokenT< edm::SimVertexContainer > simVertexToken
bool isHitOnLayer(const TrajectorySeedHitCandidate &trackerRecHit, const TrackingLayer &layer) const
tuple conf
Definition: dbtoconf.py:185
bool pass2HitsCuts(const TrajectorySeedHitCandidate &hit1, const TrajectorySeedHitCandidate &hit2) const
edm::EDGetTokenT< edm::SimTrackContainer > simTrackToken
virtual bool passSimTrackQualityCuts(const SimTrack &theSimTrack, const SimVertex &theSimVertex) const
method checks if a SimTrack fulfills the quality requirements.
std::vector< std::vector< TrackingLayer > > seedingLayers
edm::EDGetTokenT< reco::VertexCollection > recoVertexToken
const TrackerGeometry * trackerGeometry
const SeedingNode< TrackingLayer > * insertHit(const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > &hitIndicesInTree, const SeedingNode< TrackingLayer > *node, unsigned int trackerHit) const
method inserts hit into the tree structure at an empty position.
std::shared_ptr< PropagatorWithMaterial > thePropagator
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken
unsigned int getIndex() const
Definition: SeedingTree.h:107
const TrackingLayer & getTrackingLayer() const
SeedingTree< TrackingLayer > _seedingTree
bool compatibleWithPrimaryVertex(const GlobalPoint &gpos1, const GlobalPoint &gpos2, double error, bool forward) const
const MagneticFieldMap * magneticFieldMap
Definition: Run.h:41