CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TwoBodyDecayTrajectoryFactory.cc
Go to the documentation of this file.
9 
11 
14 
16 
19 
20 
23 
25 
32 {
33 public:
34 
38 
41 
45  const reco::BeamSpot &beamSpot) const;
46 
49  const ExternalPredictionCollection &external,
50  const reco::BeamSpot &beamSpot) const;
51 
52  virtual TwoBodyDecayTrajectoryFactory* clone() const { return new TwoBodyDecayTrajectoryFactory(*this); }
53 
54 protected:
56  const TwoBodyDecay &tbd,
57  const MagneticField *magField,
58  const reco::BeamSpot &beamSpot,
59  bool setParameterErrors) const;
60 
63 
65 
67 
69 
73 
74 };
75 
79 
81  : TrajectoryFactoryBase( config ),
82  theFitter( new TwoBodyDecayFitter( config ) ),
83  theNSigmaCutValue( config.getParameter< double >( "NSigmaCut" ) ),
84  theUseRefittedStateFlag( config.getParameter< bool >( "UseRefittedState" ) ),
85  theConstructTsosWithErrorsFlag( config.getParameter< bool >( "ConstructTsosWithErrors" ) )
86 {
87  produceVirtualMeasurement( config );
88 }
89 
90 
94  const reco::BeamSpot &beamSpot) const
95 {
97 
98  edm::ESHandle< MagneticField > magneticField;
99  setup.get< IdealMagneticFieldRecord >().get( magneticField );
100 
101  if ( tracks.size() == 2 )
102  {
103  // produce transient tracks from persistent tracks
104  std::vector< reco::TransientTrack > transientTracks( 2 );
105 
106  transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
107  transientTracks[0].setES( setup );
108 
109  transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
110  transientTracks[1].setES( setup );
111 
112  // estimate the decay parameters
113  TwoBodyDecay tbd = theFitter->estimate( transientTracks, theVM );
114 
115  if ( !tbd.isValid() )
116  {
117  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
118  return trajectories;
119  }
120 
121  return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, false );
122  }
123  else
124  {
125  edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
126  << "Need 2 tracks, got " << tracks.size() << ".\n";
127  }
128 
129  return trajectories;
130 }
131 
132 
136  const ExternalPredictionCollection &external,
137  const reco::BeamSpot &beamSpot) const
138 {
140 
141  edm::ESHandle< MagneticField > magneticField;
142  setup.get< IdealMagneticFieldRecord >().get( magneticField );
143 
144  if ( tracks.size() == 2 && external.size() == 2 )
145  {
146  if ( external[0].isValid() && external[1].isValid() ) // Include external estimates
147  {
148  // produce transient tracks from persistent tracks
149  std::vector< reco::TransientTrack > transientTracks( 2 );
150 
151  transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
152  transientTracks[0].setES( setup );
153 
154  transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
155  transientTracks[1].setES( setup );
156 
157  // estimate the decay parameters. the transient tracks are not really associated to the
158  // the external tsos, but this is o.k., because the only information retrieved from them
159  // is the magnetic field.
160  TwoBodyDecay tbd = theFitter->estimate( transientTracks, external, theVM );
161 
162  if ( !tbd.isValid() )
163  {
164  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
165  return trajectories;
166  }
167 
168  return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, true );
169  }
170  else
171  {
172  // Return without external estimate
173  trajectories = this->trajectories(setup, tracks, beamSpot);
174  }
175  }
176  else
177  {
178  edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
179  << "Need 2 tracks, got " << tracks.size() << ".\n";
180  }
181 
182  return trajectories;
183 }
184 
185 
188  const TwoBodyDecay& tbd,
189  const MagneticField* magField,
190  const reco::BeamSpot &beamSpot,
191  bool setParameterErrors ) const
192 {
194 
195  // get innermost valid trajectory state and hits from the tracks
196  TrajectoryInput input1 = this->innermostStateAndRecHits( tracks[0] );
197  TrajectoryInput input2 = this->innermostStateAndRecHits( tracks[1] );
198 
199  if ( !( input1.first.isValid() && input2.first.isValid() ) ) return trajectories;
200 
201  // produce TwoBodyDecayTrajectoryState (input for TwoBodyDecayTrajectory)
202  TsosContainer tsos( input1.first, input2.first );
203  ConstRecHitCollection recHits( input1.second, input2.second );
204  TwoBodyDecayTrajectoryState trajectoryState( tsos, tbd, theVM.secondaryMass(), magField );
205 
206  if ( !trajectoryState.isValid() )
207  {
208  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
209  return trajectories;
210  }
211 
212  // always use the refitted trajectory state for matching
213  TransientTrackingRecHit::ConstRecHitPointer updatedRecHit1( recHits.first.front()->clone( tsos.first ) );
214  bool valid1 = match( trajectoryState.trajectoryStates( true ).first,
215  updatedRecHit1 );
216 
217  TransientTrackingRecHit::ConstRecHitPointer updatedRecHit2( recHits.second.front()->clone( tsos.second ) );
218  bool valid2 = match( trajectoryState.trajectoryStates( true ).second,
219  updatedRecHit2 );
220 
221  if ( !valid1 || !valid2 )
222  {
223  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
224  return trajectories;
225  }
226 
227  // set the flag for reversing the RecHits to false, since they are already in the correct order.
228  TwoBodyDecayTrajectory* result = new TwoBodyDecayTrajectory( trajectoryState, recHits, magField,
230  false, beamSpot, theUseRefittedStateFlag,
232  if ( setParameterErrors && tbd.hasError() ) result->setParameterErrors( tbd.covariance() );
233  trajectories.push_back( ReferenceTrajectoryPtr( result ) );
234  return trajectories;
235 }
236 
237 
239  const TransientTrackingRecHit::ConstRecHitPointer& recHit ) const
240 {
241  LocalPoint lp1 = state.localPosition();
242  LocalPoint lp2 = recHit->localPosition();
243 
244  double deltaX = lp1.x() - lp2.x();
245  double deltaY = lp1.y() - lp2.y();
246 
247  LocalError le = recHit->localPositionError();
248 
249  double varX = le.xx();
250  double varY = le.yy();
251 
252  AlignmentPositionError* gape = recHit->det()->alignmentPositionError();
253  if ( gape )
254  {
256  LocalError lape = eft.transform( gape->globalError(), recHit->det()->surface() );
257 
258  varX += lape.xx();
259  varY += lape.yy();
260  }
261 
262  return ( ( fabs(deltaX)/sqrt(varX) < theNSigmaCutValue ) && ( fabs(deltaY)/sqrt(varY) < theNSigmaCutValue ) );
263 }
264 
265 
267 {
268  const edm::ParameterSet bsc = config.getParameter< edm::ParameterSet >( "BeamSpot" );
269  const edm::ParameterSet ppc = config.getParameter< edm::ParameterSet >( "ParticleProperties" );
270  // FIXME: Should get 3D beamspot and errors from BeamSpot input from
271  // event with extrapolation of tracks to beam line?
272 
273  GlobalPoint theBeamSpot( bsc.getParameter< double >( "MeanX" ),
274  bsc.getParameter< double >( "MeanY" ),
275  bsc.getParameter< double >( "MeanZ" ) );
276 
277  GlobalError theBeamSpotError( bsc.getParameter< double >( "VarXX" ),
278  bsc.getParameter< double >( "VarXY" ),
279  bsc.getParameter< double >( "VarYY" ),
280  bsc.getParameter< double >( "VarXZ" ),
281  bsc.getParameter< double >( "VarYZ" ),
282  bsc.getParameter< double >( "VarZZ" ) );
283 
284  theVM = VirtualMeasurement( ppc.getParameter< double >( "PrimaryMass" ),
285  ppc.getParameter< double >( "PrimaryWidth" ),
286  ppc.getParameter< double >( "SecondaryMass" ),
287  theBeamSpot, theBeamSpotError );
288  return;
289 }
290 
291 
292 DEFINE_EDM_PLUGIN( TrajectoryFactoryPlugin, TwoBodyDecayTrajectoryFactory, "TwoBodyDecayTrajectoryFactory" );
T getParameter(std::string const &) const
float xx() const
Definition: LocalError.h:19
GlobalError transform(const LocalError &le, const Surface &surf) const
TwoBodyDecayTrajectory::ConstRecHitCollection ConstRecHitCollection
TwoBodyDecayVirtualMeasurement VirtualMeasurement
void setParameterErrors(const AlgebraicSymMatrix &error)
const ReferenceTrajectoryCollection constructTrajectories(const ConstTrajTrackPairCollection &tracks, const TwoBodyDecay &tbd, const MagneticField *magField, const reco::BeamSpot &beamSpot, bool setParameterErrors) const
virtual const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const
Produce the trajectories.
T y() const
Definition: PV3DBase.h:57
std::pair< ConstRecHitContainer, ConstRecHitContainer > ConstRecHitCollection
#define input2
Definition: AMPTWrapper.h:149
bool match(const TrajectoryStateOnSurface &state, const TransientTrackingRecHit::ConstRecHitPointer &recHit) const
virtual const TwoBodyDecay estimate(const std::vector< reco::TransientTrack > &tracks, const TwoBodyDecayVirtualMeasurement &vm) const
U second(std::pair< T, U > const &p)
const AlgebraicSymMatrix & covariance(void) const
Definition: TwoBodyDecay.h:33
const bool isValid(void) const
Definition: TwoBodyDecay.h:45
float yy() const
Definition: LocalError.h:21
T sqrt(T t)
Definition: SSEVec.h:28
const bool hasError(void) const
Definition: TwoBodyDecay.h:41
tuple result
Definition: query.py:137
void produceVirtualMeasurement(const edm::ParameterSet &config)
std::pair< TrajectoryStateOnSurface, TrajectoryStateOnSurface > TsosContainer
#define input1
Definition: AMPTWrapper.h:129
const MaterialEffects materialEffects(void) const
void setES(const edm::EventSetup &es)
virtual const TrajectoryInput innermostStateAndRecHits(const ConstTrajTrackPair &track) const
AlignmentAlgorithmBase::ConstTrajTrackPairCollection ConstTrajTrackPairCollection
std::vector< ReferenceTrajectoryPtr > ReferenceTrajectoryCollection
tuple tracks
Definition: testEve_cfg.py:39
virtual TwoBodyDecayTrajectoryFactory * clone() const
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
const GlobalError & globalError() const
char state
Definition: procUtils.cc:75
TwoBodyDecayTrajectoryFactory(const edm::ParameterSet &config)
tuple config
Definition: cmsDriver.py:17
const PropagationDirection propagationDirection(void) const
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::pair< TrajectoryStateOnSurface, TransientTrackingRecHit::ConstRecHitContainer > TrajectoryInput
std::vector< TrajectoryStateOnSurface > ExternalPredictionCollection
T x() const
Definition: PV3DBase.h:56
ReferenceTrajectoryBase::ReferenceTrajectoryPtr ReferenceTrajectoryPtr
TwoBodyDecayTrajectoryState::TsosContainer TsosContainer