CMS 3D CMS Logo

TwoBodyDecayTrajectoryFactory.cc
Go to the documentation of this file.
9 
11 
13 
15 
18 
19 
22 
24 
31 {
32 public:
33 
37 
40 
44  const reco::BeamSpot &beamSpot) const override;
45 
49  const reco::BeamSpot &beamSpot) const override;
50 
51  TwoBodyDecayTrajectoryFactory* clone() const override { return new TwoBodyDecayTrajectoryFactory(*this); }
52 
53 protected:
55  const TwoBodyDecay &tbd,
56  const MagneticField *magField,
57  const reco::BeamSpot &beamSpot,
58  bool setParameterErrors) const;
59 
60  bool match(const TrajectoryStateOnSurface &state,
62 
64 
68 
73 
74 };
75 
79 
81  : TrajectoryFactoryBase(config, 2),
82  theFitter(config)
83 {
84  const edm::ParameterSet ppc = config.getParameter< edm::ParameterSet >( "ParticleProperties" );
85  thePrimaryMass = ppc.getParameter< double >( "PrimaryMass" );
86  thePrimaryWidth = ppc.getParameter< double >( "PrimaryWidth" );
87  theSecondaryMass = ppc.getParameter< double >( "SecondaryMass" );
88 
89  theNSigmaCutValue = config.getParameter< double >( "NSigmaCut" );
90  theChi2CutValue = config.getParameter< double >( "Chi2Cut" );
91  theUseRefittedStateFlag = config.getParameter< bool >( "UseRefittedState" );
92  theConstructTsosWithErrorsFlag = config.getParameter< bool >( "ConstructTsosWithErrors" );
93 }
94 
96 {
97 }
98 
102  const reco::BeamSpot &beamSpot) const
103 {
105 
107  setup.get< IdealMagneticFieldRecord >().get( magneticField );
108 
109  if ( tracks.size() == 2 )
110  {
111  // produce transient tracks from persistent tracks
112  std::vector< reco::TransientTrack > transientTracks( 2 );
113 
114  transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
115  transientTracks[0].setES( setup );
116 
117  transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
118  transientTracks[1].setES( setup );
119 
120  // estimate the decay parameters
122  TwoBodyDecay tbd = theFitter.estimate( transientTracks, vm );
123 
124  if ( !tbd.isValid() || ( tbd.chi2() > theChi2CutValue ) )
125  {
126  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
127  return trajectories;
128  }
129 
130  return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, false );
131  }
132  else
133  {
134  edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
135  << "Need 2 tracks, got " << tracks.size() << ".\n";
136  }
137 
138  return trajectories;
139 }
140 
141 
146  const reco::BeamSpot &beamSpot) const
147 {
149 
151  setup.get< IdealMagneticFieldRecord >().get( magneticField );
152 
153  if ( tracks.size() == 2 && external.size() == 2 )
154  {
155  if ( external[0].isValid() && external[1].isValid() ) // Include external estimates
156  {
157  // produce transient tracks from persistent tracks
158  std::vector< reco::TransientTrack > transientTracks( 2 );
159 
160  transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
161  transientTracks[0].setES( setup );
162 
163  transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
164  transientTracks[1].setES( setup );
165 
166  // estimate the decay parameters. the transient tracks are not really associated to the
167  // the external tsos, but this is o.k., because the only information retrieved from them
168  // is the magnetic field.
170  TwoBodyDecay tbd = theFitter.estimate(transientTracks, external, vm);
171 
172  if ( !tbd.isValid() || ( tbd.chi2() > theChi2CutValue ) )
173  {
174  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
175  return trajectories;
176  }
177 
178  return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, true );
179  }
180  else
181  {
182  // Return without external estimate
183  trajectories = this->trajectories(setup, tracks, beamSpot);
184  }
185  }
186  else
187  {
188  edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
189  << "Need 2 tracks, got " << tracks.size() << ".\n";
190  }
191 
192  return trajectories;
193 }
194 
195 
198  const TwoBodyDecay& tbd,
199  const MagneticField* magField,
200  const reco::BeamSpot &beamSpot,
201  bool setParameterErrors ) const
202 {
204 
205  // get innermost valid trajectory state and hits from the tracks
206  TrajectoryInput input1 = this->innermostStateAndRecHits( tracks[0] );
207  TrajectoryInput input2 = this->innermostStateAndRecHits( tracks[1] );
208 
209  if ( !( input1.first.isValid() && input2.first.isValid() ) ) return trajectories;
210 
211  // produce TwoBodyDecayTrajectoryState (input for TwoBodyDecayTrajectory)
212  TsosContainer tsos( input1.first, input2.first );
213  ConstRecHitCollection recHits( input1.second, input2.second );
214  TwoBodyDecayTrajectoryState trajectoryState( tsos, tbd, theSecondaryMass, magField );
215 
216  if ( !trajectoryState.isValid() )
217  {
218  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
219  return trajectories;
220  }
221 
222  // always use the refitted trajectory state for matching
223  // FIXME FIXME CLONE
224  //TrackingRecHit::ConstRecHitPointer updatedRecHit1( recHits.first.front()->clone( tsos.first ) );
225  // TrackingRecHit::ConstRecHitPointer updatedRecHit2( recHits.second.front()->clone( tsos.second ) );
226 
227  TrackingRecHit::ConstRecHitPointer updatedRecHit1( recHits.first.front() );
228  TrackingRecHit::ConstRecHitPointer updatedRecHit2( recHits.second.front() );
229 
230 
231  bool valid1 = match( trajectoryState.trajectoryStates( true ).first,
232  updatedRecHit1 );
233 
234  bool valid2 = match( trajectoryState.trajectoryStates( true ).second,
235  updatedRecHit2 );
236 
237  if ( !valid1 || !valid2 )
238  {
239  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
240  return trajectories;
241  }
242 
244  config.useBeamSpot = useBeamSpot_;
245  config.includeAPEs = includeAPEs_;
249  // set the flag for reversing the RecHits to false, since they are already in the correct order.
250  config.hitsAreReverse = false;
252  recHits,
253  magField,
254  beamSpot,
255  config);
256  if ( setParameterErrors && tbd.hasError() ) result->setParameterErrors( tbd.covariance() );
257  trajectories.push_back( ReferenceTrajectoryPtr( result ) );
258  return trajectories;
259 }
260 
261 
264 {
265  LocalPoint lp1 = state.localPosition();
266  LocalPoint lp2 = recHit->localPosition();
267 
268  double deltaX = lp1.x() - lp2.x();
269  double deltaY = lp1.y() - lp2.y();
270 
271  LocalError le = recHit->localPositionError();
272 
273  double varX = le.xx();
274  double varY = le.yy();
275 
276  return ( ( fabs(deltaX)/sqrt(varX) < theNSigmaCutValue ) && ( fabs(deltaY)/sqrt(varY) < theNSigmaCutValue ) );
277 }
278 
279 
280 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:287
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
const T & get() const
Definition: EventSetup.h:55
double chi2(void) const
Definition: TwoBodyDecay.h:47
TwoBodyDecayTrajectoryFactory(const edm::ParameterSet &config)
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