CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DualBzeroTrajectoryFactory.cc
Go to the documentation of this file.
8 
9 #include <algorithm>
10 
12 
14 
16 
17 
19 {
20 public:
23 
27  const reco::BeamSpot &beamSpot) const override;
28 
30  const ConstTrajTrackPairCollection &tracks,
31  const ExternalPredictionCollection &external,
32  const reco::BeamSpot &beamSpot) const override;
33 
34  virtual DualBzeroTrajectoryFactory* clone() const override { return new DualBzeroTrajectoryFactory(*this); }
35 
36 protected:
38  {
42  };
43 
45 
47  const Surface &surface,
48  const MagneticField *magField) const;
49 
50  double theMass;
52 
53 };
54 
58 
60  TrajectoryFactoryBase( config )
61 {
62  theMass = config.getParameter< double >( "ParticleMass" );
63  theMomentumEstimate = config.getParameter< double >( "MomentumEstimate" );
64 }
65 
66 
68 
69 
73  const reco::BeamSpot &beamSpot) const
74 {
76 
78  setup.get< IdealMagneticFieldRecord >().get( magneticField );
79 
80  ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
81 
82  while ( itTracks != tracks.end() )
83  {
84  const DualBzeroTrajectoryInput input = this->referenceStateAndRecHits( *itTracks );
85  // Check input: If all hits were rejected, the TSOS is initialized as invalid.
86  if ( input.refTsos.isValid() )
87  {
89  input.fwdRecHits,
90  input.bwdRecHits,
91  magneticField.product(),
94  theMass,
97  trajectories.push_back( ptr );
98  }
99 
100  ++itTracks;
101  }
102 
103  return trajectories;
104 }
105 
109  const ExternalPredictionCollection &external,
110  const reco::BeamSpot &beamSpot) const
111 {
113 
114  if ( tracks.size() != external.size() )
115  {
116  edm::LogInfo("ReferenceTrajectories") << "@SUB=DualBzeroTrajectoryFactory::trajectories"
117  << "Inconsistent input:\n"
118  << "\tnumber of tracks = " << tracks.size()
119  << "\tnumber of external predictions = " << external.size();
120  return trajectories;
121  }
122 
124  setup.get< IdealMagneticFieldRecord >().get( magneticField );
125 
126  ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
127  ExternalPredictionCollection::const_iterator itExternal = external.begin();
128 
129  while ( itTracks != tracks.end() )
130  {
132  // Check input: If all hits were rejected, the TSOS is initialized as invalid.
133  if ( input.refTsos.isValid() )
134  {
135  if ( (*itExternal).isValid() )
136  {
137  TrajectoryStateOnSurface propExternal =
138  propagateExternal( *itExternal, input.refTsos.surface(), magneticField.product() );
139 
140  if ( !propExternal.isValid() ) continue;
141 
142  // set the flag for reversing the RecHits to false, since they are already in the correct order.
144  input.fwdRecHits,
145  input.bwdRecHits,
146  magneticField.product(),
147  materialEffects(),
149  theMass,
152 
153  AlgebraicSymMatrix externalParamErrors( asHepMatrix<5>( propExternal.localError().matrix() ) );
154  ptr->setParameterErrors( externalParamErrors.sub( 2, 5 ) );
155  trajectories.push_back( ptr );
156  }
157  else
158  {
159 // GF: Why is the following commented? That is different from the other factories
160 // that usually have the non-external prediction ReferenceTrajectory as fall back...
161 // ReferenceTrajectoryPtr ptr( new DualBzeroReferenceTrajectory( input.refTsos,
162 // input.fwdRecHits,
163 // input.bwdRecHits,
164 // magneticField.product(),
165 // materialEffects(),
166 // propagationDirection(),
167 // theMass,
168 // theMomentumEstimate,
169 // beamSpot ) );
171  input.fwdRecHits,
172  input.bwdRecHits,
173  magneticField.product(),
174  materialEffects(),
176  theMass,
179 
180  //trajectories.push_back( ptr );
181  }
182  }
183 
184  ++itTracks;
185  ++itExternal;
186  }
187 
188  return trajectories;
189 }
190 
191 
194 {
196 
197  // get the trajectory measurements in the correct order, i.e. reverse if needed
198  Trajectory::DataContainer allTrajMeas = this->orderedTrajectoryMeasurements( *track.first );
199  Trajectory::DataContainer usedTrajMeas;
200  Trajectory::DataContainer::iterator itM;
201  // get all relevant trajectory measurements
202  for ( itM = allTrajMeas.begin(); itM != allTrajMeas.end(); itM++ )
203  {
204  if ( useRecHit( ( *itM ).recHit() ) ) usedTrajMeas.push_back( *itM );
205  }
206 
207  unsigned int iMeas = 0;
208  unsigned int nMeas = usedTrajMeas.size();
209  unsigned int nRefStateMeas = nMeas/2;
210  // get the valid RecHits
211  for ( itM = usedTrajMeas.begin(); itM != usedTrajMeas.end(); itM++, iMeas++ )
212  {
213  TransientTrackingRecHit::ConstRecHitPointer aRecHit = ( *itM ).recHit();
214 
215  if ( iMeas < nRefStateMeas ) {
216  input.bwdRecHits.push_back( aRecHit );
217  } else if ( iMeas > nRefStateMeas ) {
218  input.fwdRecHits.push_back( aRecHit );
219  } else { // iMeas == nRefStateMeas
220  if ( ( *itM ).updatedState().isValid() )
221  {
222  input.refTsos = ( *itM ).updatedState();
223  input.bwdRecHits.push_back( aRecHit );
224  input.fwdRecHits.push_back( aRecHit );
225  } else {
226  // if the tsos of the middle hit is not valid, try the next one ...
227  nRefStateMeas++;
228  input.bwdRecHits.push_back( aRecHit );
229  }
230  }
231  }
232 
233  // bring input.fwdRecHits into correct order
234  std::reverse( input.bwdRecHits.begin(), input.bwdRecHits.end() );
235 
236  return input;
237 }
238 
241  const Surface& surface,
242  const MagneticField* magField ) const
243 {
245  const std::pair< TrajectoryStateOnSurface, double > tsosWithPath =
246  propagator.propagateWithPath( external, surface );
247  return tsosWithPath.first;
248 }
249 
250 
T getParameter(std::string const &) const
const TrajectoryStateOnSurface propagateExternal(const TrajectoryStateOnSurface &external, const Surface &surface, const MagneticField *magField) const
tuple propagator
TransientTrackingRecHit::ConstRecHitContainer bwdRecHits
MaterialEffects materialEffects(void) const
const DualBzeroTrajectoryInput referenceStateAndRecHits(const ConstTrajTrackPair &track) const
virtual const Trajectory::DataContainer orderedTrajectoryMeasurements(const Trajectory &trajectory) const
virtual const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const override
Produce the reference trajectories.
AlignmentAlgorithmBase::ConstTrajTrackPair ConstTrajTrackPair
tuple magneticField
bool useRecHit(const TransientTrackingRecHit::ConstRecHitPointer &hitPtr) const
std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &fts, const Plane &plane) const override
propagation to plane with path length
static std::string const input
Definition: EdmProvDump.cc:44
const SurfaceType & surface() const
std::vector< TrajectoryMeasurement > DataContainer
Definition: Trajectory.h:44
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
const AlgebraicSymMatrix55 & matrix() const
const LocalTrajectoryError & localError() const
std::vector< ConstRecHitPointer > ConstRecHitContainer
AlignmentAlgorithmBase::ConstTrajTrackPairCollection ConstTrajTrackPairCollection
std::vector< ReferenceTrajectoryPtr > ReferenceTrajectoryCollection
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
PropagationDirection propagationDirection(void) const
CLHEP::HepSymMatrix AlgebraicSymMatrix
virtual DualBzeroTrajectoryFactory * clone() const override
A factory that produces instances of class ReferenceTrajectory from a given TrajTrackPairCollection.
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< TrajectoryStateOnSurface > ExternalPredictionCollection
TransientTrackingRecHit::ConstRecHitContainer fwdRecHits
DualBzeroTrajectoryFactory(const edm::ParameterSet &config)
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")