CMS 3D CMS Logo

TwoBodyDecayTrajectoryFactory.cc
Go to the documentation of this file.
8 
10 
12 
14 
17 
18 
21 
23 
30 {
31 public:
32 
36 
39 
43  const reco::BeamSpot &beamSpot) const override;
44 
48  const reco::BeamSpot &beamSpot) const override;
49 
50  TwoBodyDecayTrajectoryFactory* clone() const override { return new TwoBodyDecayTrajectoryFactory(*this); }
51 
52 protected:
54  const TwoBodyDecay &tbd,
55  const MagneticField *magField,
56  const reco::BeamSpot &beamSpot,
57  bool setParameterErrors) const;
58 
59  bool match(const TrajectoryStateOnSurface &state,
61 
63 
67 
72 
73 };
74 
78 
80  : TrajectoryFactoryBase(config, 2),
81  theFitter(config)
82 {
83  const edm::ParameterSet ppc = config.getParameter< edm::ParameterSet >( "ParticleProperties" );
84  thePrimaryMass = ppc.getParameter< double >( "PrimaryMass" );
85  thePrimaryWidth = ppc.getParameter< double >( "PrimaryWidth" );
86  theSecondaryMass = ppc.getParameter< double >( "SecondaryMass" );
87 
88  theNSigmaCutValue = config.getParameter< double >( "NSigmaCut" );
89  theChi2CutValue = config.getParameter< double >( "Chi2Cut" );
90  theUseRefittedStateFlag = config.getParameter< bool >( "UseRefittedState" );
91  theConstructTsosWithErrorsFlag = config.getParameter< bool >( "ConstructTsosWithErrors" );
92 }
93 
95 {
96 }
97 
101  const reco::BeamSpot &beamSpot) const
102 {
104 
106  setup.get< IdealMagneticFieldRecord >().get( magneticField );
107 
108  if ( tracks.size() == 2 )
109  {
110  // produce transient tracks from persistent tracks
111  std::vector< reco::TransientTrack > transientTracks( 2 );
112 
113  transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
114  transientTracks[0].setES( setup );
115 
116  transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
117  transientTracks[1].setES( setup );
118 
119  // estimate the decay parameters
121  TwoBodyDecay tbd = theFitter.estimate( transientTracks, vm );
122 
123  if ( !tbd.isValid() || ( tbd.chi2() > theChi2CutValue ) )
124  {
125  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
126  return trajectories;
127  }
128 
129  return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, false );
130  }
131  else
132  {
133  edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
134  << "Need 2 tracks, got " << tracks.size() << ".\n";
135  }
136 
137  return trajectories;
138 }
139 
140 
145  const reco::BeamSpot &beamSpot) const
146 {
148 
150  setup.get< IdealMagneticFieldRecord >().get( magneticField );
151 
152  if ( tracks.size() == 2 && external.size() == 2 )
153  {
154  if ( external[0].isValid() && external[1].isValid() ) // Include external estimates
155  {
156  // produce transient tracks from persistent tracks
157  std::vector< reco::TransientTrack > transientTracks( 2 );
158 
159  transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
160  transientTracks[0].setES( setup );
161 
162  transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
163  transientTracks[1].setES( setup );
164 
165  // estimate the decay parameters. the transient tracks are not really associated to the
166  // the external tsos, but this is o.k., because the only information retrieved from them
167  // is the magnetic field.
169  TwoBodyDecay tbd = theFitter.estimate(transientTracks, external, vm);
170 
171  if ( !tbd.isValid() || ( tbd.chi2() > theChi2CutValue ) )
172  {
173  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
174  return trajectories;
175  }
176 
177  return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, true );
178  }
179  else
180  {
181  // Return without external estimate
182  trajectories = this->trajectories(setup, tracks, beamSpot);
183  }
184  }
185  else
186  {
187  edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
188  << "Need 2 tracks, got " << tracks.size() << ".\n";
189  }
190 
191  return trajectories;
192 }
193 
194 
197  const TwoBodyDecay& tbd,
198  const MagneticField* magField,
199  const reco::BeamSpot &beamSpot,
200  bool setParameterErrors ) const
201 {
203 
204  // get innermost valid trajectory state and hits from the tracks
205  TrajectoryInput input1 = this->innermostStateAndRecHits( tracks[0] );
206  TrajectoryInput input2 = this->innermostStateAndRecHits( tracks[1] );
207 
208  if ( !( input1.first.isValid() && input2.first.isValid() ) ) return trajectories;
209 
210  // produce TwoBodyDecayTrajectoryState (input for TwoBodyDecayTrajectory)
211  TsosContainer tsos( input1.first, input2.first );
212  ConstRecHitCollection recHits( input1.second, input2.second );
213  TwoBodyDecayTrajectoryState trajectoryState( tsos, tbd, theSecondaryMass, magField );
214 
215  if ( !trajectoryState.isValid() )
216  {
217  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
218  return trajectories;
219  }
220 
221  // always use the refitted trajectory state for matching
222  // FIXME FIXME CLONE
223  //TrackingRecHit::ConstRecHitPointer updatedRecHit1( recHits.first.front()->clone( tsos.first ) );
224  // TrackingRecHit::ConstRecHitPointer updatedRecHit2( recHits.second.front()->clone( tsos.second ) );
225 
226  TrackingRecHit::ConstRecHitPointer updatedRecHit1( recHits.first.front() );
227  TrackingRecHit::ConstRecHitPointer updatedRecHit2( recHits.second.front() );
228 
229 
230  bool valid1 = match( trajectoryState.trajectoryStates( true ).first,
231  updatedRecHit1 );
232 
233  bool valid2 = match( trajectoryState.trajectoryStates( true ).second,
234  updatedRecHit2 );
235 
236  if ( !valid1 || !valid2 )
237  {
238  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
239  return trajectories;
240  }
241 
243  config.useBeamSpot = useBeamSpot_;
244  config.includeAPEs = includeAPEs_;
248  // set the flag for reversing the RecHits to false, since they are already in the correct order.
249  config.hitsAreReverse = false;
251  recHits,
252  magField,
253  beamSpot,
254  config);
255  if ( setParameterErrors && tbd.hasError() ) result->setParameterErrors( tbd.covariance() );
256  trajectories.push_back( ReferenceTrajectoryPtr( result ) );
257  return trajectories;
258 }
259 
260 
263 {
264  LocalPoint lp1 = state.localPosition();
265  LocalPoint lp2 = recHit->localPosition();
266 
267  double deltaX = lp1.x() - lp2.x();
268  double deltaY = lp1.y() - lp2.y();
269 
270  LocalError le = recHit->localPositionError();
271 
272  double varX = le.xx();
273  double varY = le.yy();
274 
275  return ( ( fabs(deltaX)/sqrt(varX) < theNSigmaCutValue ) && ( fabs(deltaY)/sqrt(varY) < theNSigmaCutValue ) );
276 }
277 
278 
279 DEFINE_EDM_PLUGIN( TrajectoryFactoryPlugin, TwoBodyDecayTrajectoryFactory, "TwoBodyDecayTrajectoryFactory" );
T getParameter(std::string const &) const
float xx() const
Definition: LocalError.h:24
bool hasError(void) const
Definition: TwoBodyDecay.h:45
const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const override
Produce the trajectories.
MaterialEffects materialEffects(void) 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
T y() const
Definition: PV3DBase.h:63
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: config.py:1
std::pair< ConstRecHitContainer, ConstRecHitContainer > ConstRecHitCollection
const TsosContainer & trajectoryStates(bool useRefittedState=true) const
#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)
config
Definition: looper.py:288
const AlgebraicSymMatrix & covariance(void) const
Definition: TwoBodyDecay.h:37
float yy() const
Definition: LocalError.h:26
TwoBodyDecayTrajectoryFactory * clone() const override
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
T sqrt(T t)
Definition: SSEVec.h:18
std::pair< TrajectoryStateOnSurface, TrajectoryStateOnSurface > TsosContainer
#define input1
Definition: AMPTWrapper.h:129
void setES(const edm::EventSetup &es)
bool isValid(void) const
Definition: TwoBodyDecay.h:49
virtual const TrajectoryInput innermostStateAndRecHits(const ConstTrajTrackPair &track) const
AlignmentAlgorithmBase::ConstTrajTrackPairCollection ConstTrajTrackPairCollection
std::vector< ReferenceTrajectoryPtr > ReferenceTrajectoryCollection
double chi2(void) const
Definition: TwoBodyDecay.h:47
TwoBodyDecayTrajectoryFactory(const edm::ParameterSet &config)
T get() const
Definition: EventSetup.h:63
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:62
ReferenceTrajectoryBase::ReferenceTrajectoryPtr ReferenceTrajectoryPtr
T const * product() const
Definition: ESHandle.h:86
TwoBodyDecayTrajectoryState::TsosContainer TsosContainer