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 
69 
74 
75 };
76 
80 
82  : TrajectoryFactoryBase( config ),
83  theFitter( new TwoBodyDecayFitter( config ) )
84 {
85  const edm::ParameterSet ppc = config.getParameter< edm::ParameterSet >( "ParticleProperties" );
86  thePrimaryMass = ppc.getParameter< double >( "PrimaryMass" );
87  thePrimaryWidth = ppc.getParameter< double >( "PrimaryWidth" );
88  theSecondaryMass = ppc.getParameter< double >( "SecondaryMass" );
89 
90  theNSigmaCutValue = config.getParameter< double >( "NSigmaCut" );
91  theChi2CutValue = config.getParameter< double >( "Chi2Cut" );
92  theUseRefittedStateFlag = config.getParameter< bool >( "UseRefittedState" );
93  theConstructTsosWithErrorsFlag = config.getParameter< bool >( "ConstructTsosWithErrors" );
94 }
95 
97 {
98  delete theFitter;
99 }
100 
104  const reco::BeamSpot &beamSpot) const
105 {
107 
108  edm::ESHandle< MagneticField > magneticField;
109  setup.get< IdealMagneticFieldRecord >().get( magneticField );
110 
111  if ( tracks.size() == 2 )
112  {
113  // produce transient tracks from persistent tracks
114  std::vector< reco::TransientTrack > transientTracks( 2 );
115 
116  transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
117  transientTracks[0].setES( setup );
118 
119  transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
120  transientTracks[1].setES( setup );
121 
122  // estimate the decay parameters
124  TwoBodyDecay tbd = theFitter->estimate( transientTracks, vm );
125 
126  if ( !tbd.isValid() || ( tbd.chi2() > theChi2CutValue ) )
127  {
128  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
129  return trajectories;
130  }
131 
132  return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, false );
133  }
134  else
135  {
136  edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
137  << "Need 2 tracks, got " << tracks.size() << ".\n";
138  }
139 
140  return trajectories;
141 }
142 
143 
147  const ExternalPredictionCollection &external,
148  const reco::BeamSpot &beamSpot) const
149 {
151 
152  edm::ESHandle< MagneticField > magneticField;
153  setup.get< IdealMagneticFieldRecord >().get( magneticField );
154 
155  if ( tracks.size() == 2 && external.size() == 2 )
156  {
157  if ( external[0].isValid() && external[1].isValid() ) // Include external estimates
158  {
159  // produce transient tracks from persistent tracks
160  std::vector< reco::TransientTrack > transientTracks( 2 );
161 
162  transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
163  transientTracks[0].setES( setup );
164 
165  transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
166  transientTracks[1].setES( setup );
167 
168  // estimate the decay parameters. the transient tracks are not really associated to the
169  // the external tsos, but this is o.k., because the only information retrieved from them
170  // is the magnetic field.
172  TwoBodyDecay tbd = theFitter->estimate( transientTracks, external, vm );
173 
174  if ( !tbd.isValid() || ( tbd.chi2() > theChi2CutValue ) )
175  {
176  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
177  return trajectories;
178  }
179 
180  return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, true );
181  }
182  else
183  {
184  // Return without external estimate
185  trajectories = this->trajectories(setup, tracks, beamSpot);
186  }
187  }
188  else
189  {
190  edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
191  << "Need 2 tracks, got " << tracks.size() << ".\n";
192  }
193 
194  return trajectories;
195 }
196 
197 
200  const TwoBodyDecay& tbd,
201  const MagneticField* magField,
202  const reco::BeamSpot &beamSpot,
203  bool setParameterErrors ) const
204 {
206 
207  // get innermost valid trajectory state and hits from the tracks
208  TrajectoryInput input1 = this->innermostStateAndRecHits( tracks[0] );
209  TrajectoryInput input2 = this->innermostStateAndRecHits( tracks[1] );
210 
211  if ( !( input1.first.isValid() && input2.first.isValid() ) ) return trajectories;
212 
213  // produce TwoBodyDecayTrajectoryState (input for TwoBodyDecayTrajectory)
214  TsosContainer tsos( input1.first, input2.first );
215  ConstRecHitCollection recHits( input1.second, input2.second );
216  TwoBodyDecayTrajectoryState trajectoryState( tsos, tbd, theSecondaryMass, magField );
217 
218  if ( !trajectoryState.isValid() )
219  {
220  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
221  return trajectories;
222  }
223 
224  // always use the refitted trajectory state for matching
225  TransientTrackingRecHit::ConstRecHitPointer updatedRecHit1( recHits.first.front()->clone( tsos.first ) );
226  bool valid1 = match( trajectoryState.trajectoryStates( true ).first,
227  updatedRecHit1 );
228 
229  TransientTrackingRecHit::ConstRecHitPointer updatedRecHit2( recHits.second.front()->clone( tsos.second ) );
230  bool valid2 = match( trajectoryState.trajectoryStates( true ).second,
231  updatedRecHit2 );
232 
233  if ( !valid1 || !valid2 )
234  {
235  trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
236  return trajectories;
237  }
238 
239  // set the flag for reversing the RecHits to false, since they are already in the correct order.
240  TwoBodyDecayTrajectory* result = new TwoBodyDecayTrajectory( trajectoryState, recHits, magField,
242  false, beamSpot, theUseRefittedStateFlag,
244  if ( setParameterErrors && tbd.hasError() ) result->setParameterErrors( tbd.covariance() );
245  trajectories.push_back( ReferenceTrajectoryPtr( result ) );
246  return trajectories;
247 }
248 
249 
251  const TransientTrackingRecHit::ConstRecHitPointer& recHit ) const
252 {
253  LocalPoint lp1 = state.localPosition();
254  LocalPoint lp2 = recHit->localPosition();
255 
256  double deltaX = lp1.x() - lp2.x();
257  double deltaY = lp1.y() - lp2.y();
258 
259  LocalError le = recHit->localPositionError();
260 
261  double varX = le.xx();
262  double varY = le.yy();
263 
264  AlignmentPositionError* gape = recHit->det()->alignmentPositionError();
265  if ( gape )
266  {
268  LocalError lape = eft.transform( gape->globalError(), recHit->det()->surface() );
269 
270  varX += lape.xx();
271  varY += lape.yy();
272  }
273 
274  return ( ( fabs(deltaX)/sqrt(varX) < theNSigmaCutValue ) && ( fabs(deltaY)/sqrt(varY) < theNSigmaCutValue ) );
275 }
276 
277 
278 DEFINE_EDM_PLUGIN( TrajectoryFactoryPlugin, TwoBodyDecayTrajectoryFactory, "TwoBodyDecayTrajectoryFactory" );
T getParameter(std::string const &) const
static GlobalError transform(const LocalError &le, const Surface &surf)
float xx() const
Definition: LocalError.h:24
bool hasError(void) const
Definition: TwoBodyDecay.h:45
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
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:63
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)
const AlgebraicSymMatrix & covariance(void) const
Definition: TwoBodyDecay.h:37
float yy() const
Definition: LocalError.h:26
T sqrt(T t)
Definition: SSEVec.h:48
tuple result
Definition: query.py:137
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
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
double chi2(void) const
Definition: TwoBodyDecay.h:47
const GlobalError & globalError() const
char state
Definition: procUtils.cc:75
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
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
TwoBodyDecayTrajectoryState::TsosContainer TsosContainer