00001
00002 #include "BzeroReferenceTrajectoryFactory.h"
00003 #include "Alignment/ReferenceTrajectories/interface/BzeroReferenceTrajectory.h"
00004 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryPlugin.h"
00005
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/Framework/interface/EventSetup.h"
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010
00011 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00012
00013
00014 BzeroReferenceTrajectoryFactory::BzeroReferenceTrajectoryFactory( const edm::ParameterSet & config ) :
00015 TrajectoryFactoryBase( config )
00016 {
00017 theMass = config.getParameter< double >( "ParticleMass" );
00018 theMomentumEstimate = config.getParameter< double >( "MomentumEstimate" );
00019 }
00020
00021
00022 BzeroReferenceTrajectoryFactory::~BzeroReferenceTrajectoryFactory( void ) {}
00023
00024
00025 const BzeroReferenceTrajectoryFactory::ReferenceTrajectoryCollection
00026 BzeroReferenceTrajectoryFactory::trajectories( const edm::EventSetup & setup,
00027 const ConstTrajTrackPairCollection & tracks ) const
00028 {
00029 ReferenceTrajectoryCollection trajectories;
00030
00031 edm::ESHandle< MagneticField > magneticField;
00032 setup.get< IdealMagneticFieldRecord >().get( magneticField );
00033
00034 ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
00035
00036 while ( itTracks != tracks.end() )
00037 {
00038 TrajectoryInput input = this->innermostStateAndRecHits( *itTracks );
00039
00040 if ( input.first.isValid() )
00041 {
00042
00043 trajectories.push_back(ReferenceTrajectoryPtr(new BzeroReferenceTrajectory(input.first, input.second, false,
00044 magneticField.product(),
00045 materialEffects(),
00046 propagationDirection(),
00047 theMass, theMomentumEstimate)));
00048 }
00049
00050 ++itTracks;
00051 }
00052
00053 return trajectories;
00054 }
00055
00056
00057 const BzeroReferenceTrajectoryFactory::ReferenceTrajectoryCollection
00058 BzeroReferenceTrajectoryFactory::trajectories( const edm::EventSetup & setup,
00059 const ConstTrajTrackPairCollection& tracks,
00060 const ExternalPredictionCollection& external ) const
00061 {
00062 ReferenceTrajectoryCollection trajectories;
00063
00064 if ( tracks.size() != external.size() )
00065 {
00066 edm::LogInfo("ReferenceTrajectories") << "@SUB=BzeroReferenceTrajectoryFactory::trajectories"
00067 << "Inconsistent input:\n"
00068 << "\tnumber of tracks = " << tracks.size()
00069 << "\tnumber of external predictions = " << external.size();
00070 return trajectories;
00071 }
00072
00073 edm::ESHandle< MagneticField > magneticField;
00074 setup.get< IdealMagneticFieldRecord >().get( magneticField );
00075
00076 ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
00077 ExternalPredictionCollection::const_iterator itExternal = external.begin();
00078
00079 while ( itTracks != tracks.end() )
00080 {
00081 TrajectoryInput input = innermostStateAndRecHits( *itTracks );
00082
00083 if ( input.first.isValid() )
00084 {
00085 if ( (*itExternal).isValid() && sameSurface( (*itExternal).surface(), input.first.surface() ) )
00086 {
00087
00088 ReferenceTrajectoryPtr refTraj( new BzeroReferenceTrajectory( *itExternal, input.second, false,
00089 magneticField.product(), materialEffects(),
00090 propagationDirection(), theMass,
00091 theMomentumEstimate ) );
00092
00093 AlgebraicSymMatrix externalParamErrors( asHepMatrix<5>( (*itExternal).localError().matrix() ) );
00094 refTraj->setParameterErrors( externalParamErrors.sub( 2, 5 ) );
00095
00096 trajectories.push_back( refTraj );
00097 }
00098 else
00099 {
00100 trajectories.push_back(ReferenceTrajectoryPtr(new BzeroReferenceTrajectory(input.first, input.second, false,
00101 magneticField.product(),
00102 materialEffects(),
00103 propagationDirection(),
00104 theMass, theMomentumEstimate)));
00105 }
00106 }
00107
00108 ++itTracks;
00109 ++itExternal;
00110 }
00111
00112 return trajectories;
00113 }
00114
00115
00116
00117 DEFINE_EDM_PLUGIN( TrajectoryFactoryPlugin, BzeroReferenceTrajectoryFactory, "BzeroReferenceTrajectoryFactory" );