![]() |
![]() |
#include <Alignment/ReferenceTrajectories/plugins/TwoBodyDecayTrajectoryFactory.h>
Definition at line 10 of file TwoBodyDecayTrajectoryFactory.h.
typedef TwoBodyDecayTrajectory::ConstRecHitCollection TwoBodyDecayTrajectoryFactory::ConstRecHitCollection |
Definition at line 17 of file TwoBodyDecayTrajectoryFactory.h.
Definition at line 16 of file TwoBodyDecayTrajectoryFactory.h.
Definition at line 15 of file TwoBodyDecayTrajectoryFactory.h.
TwoBodyDecayTrajectoryFactory::TwoBodyDecayTrajectoryFactory | ( | const edm::ParameterSet & | config | ) |
Definition at line 23 of file TwoBodyDecayTrajectoryFactory.cc.
References produceVirtualMeasurement().
Referenced by clone().
00024 : TrajectoryFactoryBase( config ), 00025 theFitter( new TwoBodyDecayFitter( config ) ), 00026 theNSigmaCutValue( config.getParameter< double >( "NSigmaCut" ) ), 00027 theUseRefittedStateFlag( config.getParameter< bool >( "UseRefittedState" ) ), 00028 theConstructTsosWithErrorsFlag( config.getParameter< bool >( "ConstructTsosWithErrors" ) ) 00029 { 00030 produceVirtualMeasurement( config ); 00031 }
TwoBodyDecayTrajectoryFactory::~TwoBodyDecayTrajectoryFactory | ( | void | ) | [inline] |
virtual TwoBodyDecayTrajectoryFactory* TwoBodyDecayTrajectoryFactory::clone | ( | void | ) | const [inline, virtual] |
Implements TrajectoryFactoryBase.
Definition at line 30 of file TwoBodyDecayTrajectoryFactory.h.
References TwoBodyDecayTrajectoryFactory().
00030 { return new TwoBodyDecayTrajectoryFactory( *this ); }
const TwoBodyDecayTrajectoryFactory::ReferenceTrajectoryCollection TwoBodyDecayTrajectoryFactory::constructTrajectories | ( | const ConstTrajTrackPairCollection & | tracks, | |
const TwoBodyDecay & | tbd, | |||
const MagneticField * | magField, | |||
bool | setParameterErrors | |||
) | const [protected] |
Definition at line 128 of file TwoBodyDecayTrajectoryFactory.cc.
References TwoBodyDecay::covariance(), TwoBodyDecay::hasError(), TrajectoryFactoryBase::innermostStateAndRecHits(), match(), TrajectoryFactoryBase::materialEffects(), TrajectoryFactoryBase::propagationDirection(), HLT_VtxMuL3::result, TwoBodyDecayVirtualMeasurement::secondaryMass(), ReferenceTrajectoryBase::setParameterErrors(), theConstructTsosWithErrorsFlag, theUseRefittedStateFlag, theVM, and trajectories().
Referenced by trajectories().
00132 { 00133 ReferenceTrajectoryCollection trajectories; 00134 00135 // get innermost valid trajectory state and hits from the tracks 00136 TrajectoryInput input1 = this->innermostStateAndRecHits( tracks[0] ); 00137 TrajectoryInput input2 = this->innermostStateAndRecHits( tracks[1] ); 00138 00139 if ( !( input1.first.isValid() && input2.first.isValid() ) ) return trajectories; 00140 00141 // produce TwoBodyDecayTrajectoryState (input for TwoBodyDecayTrajectory) 00142 TsosContainer tsos( input1.first, input2.first ); 00143 ConstRecHitCollection recHits( input1.second, input2.second ); 00144 TwoBodyDecayTrajectoryState trajectoryState( tsos, tbd, theVM.secondaryMass(), magField ); 00145 00146 if ( !trajectoryState.isValid() ) 00147 { 00148 trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) ); 00149 return trajectories; 00150 } 00151 00152 // always use the refitted trajectory state for matching 00153 TransientTrackingRecHit::ConstRecHitPointer updatedRecHit1( recHits.first.front()->clone( tsos.first ) ); 00154 bool valid1 = match( trajectoryState.trajectoryStates( true ).first, 00155 updatedRecHit1 ); 00156 00157 TransientTrackingRecHit::ConstRecHitPointer updatedRecHit2( recHits.second.front()->clone( tsos.second ) ); 00158 bool valid2 = match( trajectoryState.trajectoryStates( true ).second, 00159 updatedRecHit2 ); 00160 00161 if ( !valid1 || !valid2 ) 00162 { 00163 trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) ); 00164 return trajectories; 00165 } 00166 00167 // set the flag for reversing the RecHits to false, since they are already in the correct order. 00168 TwoBodyDecayTrajectory* result = new TwoBodyDecayTrajectory( trajectoryState, recHits, magField, 00169 materialEffects(), propagationDirection(), 00170 false, theUseRefittedStateFlag, 00171 theConstructTsosWithErrorsFlag ); 00172 if ( setParameterErrors && tbd.hasError() ) result->setParameterErrors( tbd.covariance() ); 00173 trajectories.push_back( ReferenceTrajectoryPtr( result ) ); 00174 return trajectories; 00175 }
bool TwoBodyDecayTrajectoryFactory::match | ( | const TrajectoryStateOnSurface & | state, | |
const TransientTrackingRecHit::ConstRecHitPointer & | recHit | |||
) | const [protected] |
Definition at line 178 of file TwoBodyDecayTrajectoryFactory.cc.
References AlignmentPositionError::globalError(), asciidump::le, TrajectoryStateOnSurface::localPosition(), lp1, funct::sqrt(), theNSigmaCutValue, ErrorFrameTransformer::transform(), PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), and LocalError::yy().
Referenced by constructTrajectories().
00180 { 00181 LocalPoint lp1 = state.localPosition(); 00182 LocalPoint lp2 = recHit->localPosition(); 00183 00184 double deltaX = lp1.x() - lp2.x(); 00185 double deltaY = lp1.y() - lp2.y(); 00186 00187 LocalError le = recHit->localPositionError(); 00188 00189 double varX = le.xx(); 00190 double varY = le.yy(); 00191 00192 AlignmentPositionError* gape = recHit->det()->alignmentPositionError(); 00193 if ( gape ) 00194 { 00195 ErrorFrameTransformer eft; 00196 LocalError lape = eft.transform( gape->globalError(), recHit->det()->surface() ); 00197 00198 varX += lape.xx(); 00199 varY += lape.yy(); 00200 } 00201 00202 return ( ( fabs(deltaX)/sqrt(varX) < theNSigmaCutValue ) && ( fabs(deltaY)/sqrt(varY) < theNSigmaCutValue ) ); 00203 }
void TwoBodyDecayTrajectoryFactory::produceVirtualMeasurement | ( | const edm::ParameterSet & | config | ) | [protected] |
Definition at line 206 of file TwoBodyDecayTrajectoryFactory.cc.
References edm::ParameterSet::getParameter(), theBeamSpot, theBeamSpotError, and theVM.
Referenced by TwoBodyDecayTrajectoryFactory().
00207 { 00208 const edm::ParameterSet bsc = config.getParameter< edm::ParameterSet >( "BeamSpot" ); 00209 const edm::ParameterSet ppc = config.getParameter< edm::ParameterSet >( "ParticleProperties" ); 00210 00211 theBeamSpot = GlobalPoint( bsc.getParameter< double >( "MeanX" ), 00212 bsc.getParameter< double >( "MeanY" ), 00213 bsc.getParameter< double >( "MeanZ" ) ); 00214 00215 theBeamSpotError = GlobalError( bsc.getParameter< double >( "VarXX" ), 00216 bsc.getParameter< double >( "VarXY" ), 00217 bsc.getParameter< double >( "VarYY" ), 00218 bsc.getParameter< double >( "VarXZ" ), 00219 bsc.getParameter< double >( "VarYZ" ), 00220 bsc.getParameter< double >( "VarZZ" ) ); 00221 00222 theVM = VirtualMeasurement( ppc.getParameter< double >( "PrimaryMass" ), 00223 ppc.getParameter< double >( "PrimaryWidth" ), 00224 ppc.getParameter< double >( "SecondaryMass" ), 00225 theBeamSpot, theBeamSpotError ); 00226 return; 00227 }
const TrajectoryFactoryBase::ReferenceTrajectoryCollection TwoBodyDecayTrajectoryFactory::trajectories | ( | const edm::EventSetup & | setup, | |
const ConstTrajTrackPairCollection & | tracks, | |||
const ExternalPredictionCollection & | external | |||
) | const [virtual] |
Definition at line 76 of file TwoBodyDecayTrajectoryFactory.cc.
References constructTrajectories(), TwoBodyDecayFitter::estimate(), edm::EventSetup::get(), TwoBodyDecay::isValid(), edm::ESHandle< T >::product(), edm::second(), theFitter, theVM, trajectories(), and true.
00079 { 00080 ReferenceTrajectoryCollection trajectories; 00081 00082 edm::ESHandle< MagneticField > magneticField; 00083 setup.get< IdealMagneticFieldRecord >().get( magneticField ); 00084 00085 if ( tracks.size() == 2 && external.size() == 2 ) 00086 { 00087 if ( external[0].isValid() && external[1].isValid() ) // Include external estimates 00088 { 00089 // produce transient tracks from persistent tracks 00090 std::vector< reco::TransientTrack > transientTracks( 2 ); 00091 00092 transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() ); 00093 transientTracks[0].setES( setup ); 00094 00095 transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() ); 00096 transientTracks[1].setES( setup ); 00097 00098 // estimate the decay parameters. the transient tracks are not really associated to the 00099 // the external tsos, but this is o.k., because the only information retrieved from them 00100 // is the magnetic field. 00101 TwoBodyDecay tbd = theFitter->estimate( transientTracks, external, theVM ); 00102 00103 if ( !tbd.isValid() ) 00104 { 00105 trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) ); 00106 return trajectories; 00107 } 00108 00109 return constructTrajectories( tracks, tbd, magneticField.product(), true ); 00110 } 00111 else 00112 { 00113 // Return without external estimate 00114 trajectories = this->trajectories( setup, tracks ); 00115 } 00116 } 00117 else 00118 { 00119 edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories" 00120 << "Need 2 tracks, got " << tracks.size() << ".\n"; 00121 } 00122 00123 return trajectories; 00124 }
const TrajectoryFactoryBase::ReferenceTrajectoryCollection TwoBodyDecayTrajectoryFactory::trajectories | ( | const edm::EventSetup & | setup, | |
const ConstTrajTrackPairCollection & | tracks | |||
) | const [virtual] |
Produce the trajectories.
Definition at line 35 of file TwoBodyDecayTrajectoryFactory.cc.
References constructTrajectories(), TwoBodyDecayFitter::estimate(), false, edm::EventSetup::get(), TwoBodyDecay::isValid(), edm::ESHandle< T >::product(), edm::second(), theFitter, and theVM.
Referenced by constructTrajectories(), and trajectories().
00037 { 00038 ReferenceTrajectoryCollection trajectories; 00039 00040 edm::ESHandle< MagneticField > magneticField; 00041 setup.get< IdealMagneticFieldRecord >().get( magneticField ); 00042 00043 if ( tracks.size() == 2 ) 00044 { 00045 // produce transient tracks from persistent tracks 00046 std::vector< reco::TransientTrack > transientTracks( 2 ); 00047 00048 transientTracks[0] = reco::TransientTrack( *tracks[0].second, magneticField.product() ); 00049 transientTracks[0].setES( setup ); 00050 00051 transientTracks[1] = reco::TransientTrack( *tracks[1].second, magneticField.product() ); 00052 transientTracks[1].setES( setup ); 00053 00054 // estimate the decay parameters 00055 TwoBodyDecay tbd = theFitter->estimate( transientTracks, theVM ); 00056 00057 if ( !tbd.isValid() ) 00058 { 00059 trajectories.push_back( ReferenceTrajectoryPtr( new TwoBodyDecayTrajectory() ) ); 00060 return trajectories; 00061 } 00062 00063 return constructTrajectories( tracks, tbd, magneticField.product(), false ); 00064 } 00065 else 00066 { 00067 edm::LogInfo( "ReferenceTrajectories" ) << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories" 00068 << "Need 2 tracks, got " << tracks.size() << ".\n"; 00069 } 00070 00071 return trajectories; 00072 }
Definition at line 45 of file TwoBodyDecayTrajectoryFactory.h.
Referenced by produceVirtualMeasurement().
Definition at line 46 of file TwoBodyDecayTrajectoryFactory.h.
Referenced by produceVirtualMeasurement().
Definition at line 52 of file TwoBodyDecayTrajectoryFactory.h.
Referenced by constructTrajectories().
double TwoBodyDecayTrajectoryFactory::theNSigmaCutValue [protected] |
Definition at line 51 of file TwoBodyDecayTrajectoryFactory.h.
Referenced by constructTrajectories().
Definition at line 44 of file TwoBodyDecayTrajectoryFactory.h.
Referenced by constructTrajectories(), produceVirtualMeasurement(), and trajectories().