00001 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00002 #include "DataFormats/GeometrySurface/interface/Surface.h"
00003 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
00004 #include "DataFormats/Math/interface/Vector.h"
00005 #include "DataFormats/Math/interface/Error.h"
00006 #include "TrackingTools/TrajectoryState/interface/CopyUsingClone.h"
00007 #include "RecoVertex/VertexTools/interface/PerigeeLinearizedTrackState.h"
00008 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryPlugin.h"
00009
00010 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00011
00012 #include "DataFormats/GeometryCommonDetAlgo/interface/ErrorFrameTransformer.h"
00013 #include "DataFormats/TrackingRecHit/interface/AlignmentPositionError.h"
00014
00015 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00016
00017 #include "FWCore/Framework/interface/ESHandle.h"
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019
00020
00021 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryBase.h"
00022 #include "Alignment/ReferenceTrajectories/interface/TwoBodyDecayTrajectory.h"
00023
00024 #include "Alignment/TwoBodyDecay/interface/TwoBodyDecayFitter.h"
00025
00031 class TwoBodyDecayTrajectoryFactory : public TrajectoryFactoryBase
00032 {
00033 public:
00034
00035 typedef TwoBodyDecayVirtualMeasurement VirtualMeasurement;
00036 typedef TwoBodyDecayTrajectoryState::TsosContainer TsosContainer;
00037 typedef TwoBodyDecayTrajectory::ConstRecHitCollection ConstRecHitCollection;
00038
00039 TwoBodyDecayTrajectoryFactory(const edm::ParameterSet &config);
00040 ~TwoBodyDecayTrajectoryFactory();
00041
00043 virtual const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup,
00044 const ConstTrajTrackPairCollection &tracks,
00045 const reco::BeamSpot &beamSpot) const;
00046
00047 virtual const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup,
00048 const ConstTrajTrackPairCollection &tracks,
00049 const ExternalPredictionCollection &external,
00050 const reco::BeamSpot &beamSpot) const;
00051
00052 virtual TwoBodyDecayTrajectoryFactory* clone() const { return new TwoBodyDecayTrajectoryFactory(*this); }
00053
00054 protected:
00055 const ReferenceTrajectoryCollection constructTrajectories(const ConstTrajTrackPairCollection &tracks,
00056 const TwoBodyDecay &tbd,
00057 const MagneticField *magField,
00058 const reco::BeamSpot &beamSpot,
00059 bool setParameterErrors) const;
00060
00061 bool match(const TrajectoryStateOnSurface &state,
00062 const TransientTrackingRecHit::ConstRecHitPointer &recHit) const;
00063
00064 TwoBodyDecayFitter* theFitter;
00065
00066 double thePrimaryMass;
00067 double thePrimaryWidth;
00068 double theSecondaryMass;
00069
00070 double theNSigmaCutValue;
00071 double theChi2CutValue;
00072 bool theUseRefittedStateFlag;
00073 bool theConstructTsosWithErrorsFlag;
00074
00075 };
00076
00080
00081 TwoBodyDecayTrajectoryFactory::TwoBodyDecayTrajectoryFactory( const edm::ParameterSet& config )
00082 : TrajectoryFactoryBase( config ),
00083 theFitter( new TwoBodyDecayFitter( config ) )
00084 {
00085 const edm::ParameterSet ppc = config.getParameter< edm::ParameterSet >( "ParticleProperties" );
00086 thePrimaryMass = ppc.getParameter< double >( "PrimaryMass" );
00087 thePrimaryWidth = ppc.getParameter< double >( "PrimaryWidth" );
00088 theSecondaryMass = ppc.getParameter< double >( "SecondaryMass" );
00089
00090 theNSigmaCutValue = config.getParameter< double >( "NSigmaCut" );
00091 theChi2CutValue = config.getParameter< double >( "Chi2Cut" );
00092 theUseRefittedStateFlag = config.getParameter< bool >( "UseRefittedState" );
00093 theConstructTsosWithErrorsFlag = config.getParameter< bool >( "ConstructTsosWithErrors" );
00094 }
00095
00096 TwoBodyDecayTrajectoryFactory::~TwoBodyDecayTrajectoryFactory()
00097 {
00098 delete theFitter;
00099 }
00100
00101 const TrajectoryFactoryBase::ReferenceTrajectoryCollection
00102 TwoBodyDecayTrajectoryFactory::trajectories(const edm::EventSetup &setup,
00103 const ConstTrajTrackPairCollection &tracks,
00104 const reco::BeamSpot &beamSpot) const
00105 {
00106 ReferenceTrajectoryCollection trajectories;
00107
00108 edm::ESHandle< MagneticField > magneticField;
00109 setup.get< IdealMagneticFieldRecord >().get( magneticField );
00110
00111 if ( tracks.size() == 2 )
00112 {
00113
00114 std::vector< reco::TransientTrack > transientTracks( 2 );
00115
00116 transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
00117 transientTracks[0].setES( setup );
00118
00119 transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
00120 transientTracks[1].setES( setup );
00121
00122
00123 VirtualMeasurement vm( thePrimaryMass, thePrimaryWidth, theSecondaryMass, beamSpot );
00124 TwoBodyDecay tbd = theFitter->estimate( transientTracks, vm );
00125
00126 if ( !tbd.isValid() || ( tbd.chi2() > theChi2CutValue ) )
00127 {
00128 trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
00129 return trajectories;
00130 }
00131
00132 return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, false );
00133 }
00134 else
00135 {
00136 edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
00137 << "Need 2 tracks, got " << tracks.size() << ".\n";
00138 }
00139
00140 return trajectories;
00141 }
00142
00143
00144 const TrajectoryFactoryBase::ReferenceTrajectoryCollection
00145 TwoBodyDecayTrajectoryFactory::trajectories(const edm::EventSetup &setup,
00146 const ConstTrajTrackPairCollection &tracks,
00147 const ExternalPredictionCollection &external,
00148 const reco::BeamSpot &beamSpot) const
00149 {
00150 ReferenceTrajectoryCollection trajectories;
00151
00152 edm::ESHandle< MagneticField > magneticField;
00153 setup.get< IdealMagneticFieldRecord >().get( magneticField );
00154
00155 if ( tracks.size() == 2 && external.size() == 2 )
00156 {
00157 if ( external[0].isValid() && external[1].isValid() )
00158 {
00159
00160 std::vector< reco::TransientTrack > transientTracks( 2 );
00161
00162 transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() );
00163 transientTracks[0].setES( setup );
00164
00165 transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() );
00166 transientTracks[1].setES( setup );
00167
00168
00169
00170
00171 VirtualMeasurement vm( thePrimaryMass, thePrimaryWidth, theSecondaryMass, beamSpot );
00172 TwoBodyDecay tbd = theFitter->estimate( transientTracks, external, vm );
00173
00174 if ( !tbd.isValid() || ( tbd.chi2() > theChi2CutValue ) )
00175 {
00176 trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
00177 return trajectories;
00178 }
00179
00180 return constructTrajectories( tracks, tbd, magneticField.product(), beamSpot, true );
00181 }
00182 else
00183 {
00184
00185 trajectories = this->trajectories(setup, tracks, beamSpot);
00186 }
00187 }
00188 else
00189 {
00190 edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
00191 << "Need 2 tracks, got " << tracks.size() << ".\n";
00192 }
00193
00194 return trajectories;
00195 }
00196
00197
00198 const TwoBodyDecayTrajectoryFactory::ReferenceTrajectoryCollection
00199 TwoBodyDecayTrajectoryFactory::constructTrajectories( const ConstTrajTrackPairCollection& tracks,
00200 const TwoBodyDecay& tbd,
00201 const MagneticField* magField,
00202 const reco::BeamSpot &beamSpot,
00203 bool setParameterErrors ) const
00204 {
00205 ReferenceTrajectoryCollection trajectories;
00206
00207
00208 TrajectoryInput input1 = this->innermostStateAndRecHits( tracks[0] );
00209 TrajectoryInput input2 = this->innermostStateAndRecHits( tracks[1] );
00210
00211 if ( !( input1.first.isValid() && input2.first.isValid() ) ) return trajectories;
00212
00213
00214 TsosContainer tsos( input1.first, input2.first );
00215 ConstRecHitCollection recHits( input1.second, input2.second );
00216 TwoBodyDecayTrajectoryState trajectoryState( tsos, tbd, theSecondaryMass, magField );
00217
00218 if ( !trajectoryState.isValid() )
00219 {
00220 trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
00221 return trajectories;
00222 }
00223
00224
00225 TransientTrackingRecHit::ConstRecHitPointer updatedRecHit1( recHits.first.front()->clone( tsos.first ) );
00226 bool valid1 = match( trajectoryState.trajectoryStates( true ).first,
00227 updatedRecHit1 );
00228
00229 TransientTrackingRecHit::ConstRecHitPointer updatedRecHit2( recHits.second.front()->clone( tsos.second ) );
00230 bool valid2 = match( trajectoryState.trajectoryStates( true ).second,
00231 updatedRecHit2 );
00232
00233 if ( !valid1 || !valid2 )
00234 {
00235 trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) );
00236 return trajectories;
00237 }
00238
00239
00240 TwoBodyDecayTrajectory* result = new TwoBodyDecayTrajectory( trajectoryState, recHits, magField,
00241 materialEffects(), propagationDirection(),
00242 false, beamSpot, theUseRefittedStateFlag,
00243 theConstructTsosWithErrorsFlag );
00244 if ( setParameterErrors && tbd.hasError() ) result->setParameterErrors( tbd.covariance() );
00245 trajectories.push_back( ReferenceTrajectoryPtr( result ) );
00246 return trajectories;
00247 }
00248
00249
00250 bool TwoBodyDecayTrajectoryFactory::match( const TrajectoryStateOnSurface& state,
00251 const TransientTrackingRecHit::ConstRecHitPointer& recHit ) const
00252 {
00253 LocalPoint lp1 = state.localPosition();
00254 LocalPoint lp2 = recHit->localPosition();
00255
00256 double deltaX = lp1.x() - lp2.x();
00257 double deltaY = lp1.y() - lp2.y();
00258
00259 LocalError le = recHit->localPositionError();
00260
00261 double varX = le.xx();
00262 double varY = le.yy();
00263
00264 AlignmentPositionError* gape = recHit->det()->alignmentPositionError();
00265 if ( gape )
00266 {
00267 ErrorFrameTransformer eft;
00268 LocalError lape = eft.transform( gape->globalError(), recHit->det()->surface() );
00269
00270 varX += lape.xx();
00271 varY += lape.yy();
00272 }
00273
00274 return ( ( fabs(deltaX)/sqrt(varX) < theNSigmaCutValue ) && ( fabs(deltaY)/sqrt(varY) < theNSigmaCutValue ) );
00275 }
00276
00277
00278 DEFINE_EDM_PLUGIN( TrajectoryFactoryPlugin, TwoBodyDecayTrajectoryFactory, "TwoBodyDecayTrajectoryFactory" );