00001
00002 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00003 #include "DataFormats/GeometrySurface/interface/Surface.h"
00004 #include "Alignment/ReferenceTrajectories/interface/TwoBodyDecayTrajectory.h"
00005 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
00006 #include "DataFormats/Math/interface/Error.h"
00007
00008 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00009
00010
00011 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00012
00013 TwoBodyDecayTrajectory::TwoBodyDecayTrajectory( const TwoBodyDecayTrajectoryState& trajectoryState,
00014 const ConstRecHitCollection & recHits,
00015 const MagneticField* magField,
00016 MaterialEffects materialEffects,
00017 PropagationDirection propDir,
00018 bool hitsAreReverse,
00019 const reco::BeamSpot &beamSpot,
00020 bool useRefittedState,
00021 bool constructTsosWithErrors )
00022
00023 : ReferenceTrajectoryBase(
00024 TwoBodyDecayParameters::dimension, recHits.first.size() + recHits.second.size(),
00025 (materialEffects >= breakPoints) ? 2*(recHits.first.size() + recHits.second.size())-4 : 0,
00026 (materialEffects >= breakPoints) ? 2*(recHits.first.size() + recHits.second.size())-4 : 0 )
00027 {
00028 if ( hitsAreReverse )
00029 {
00030 TransientTrackingRecHit::ConstRecHitContainer::const_reverse_iterator itRecHits;
00031 ConstRecHitCollection fwdRecHits;
00032
00033 fwdRecHits.first.reserve( recHits.first.size() );
00034 for ( itRecHits = recHits.first.rbegin(); itRecHits != recHits.first.rend(); ++itRecHits )
00035 {
00036 fwdRecHits.first.push_back( *itRecHits );
00037 }
00038
00039 fwdRecHits.second.reserve( recHits.second.size() );
00040 for ( itRecHits = recHits.second.rbegin(); itRecHits != recHits.second.rend(); ++itRecHits )
00041 {
00042 fwdRecHits.second.push_back( *itRecHits );
00043 }
00044
00045 theValidityFlag = this->construct( trajectoryState, fwdRecHits, magField, materialEffects, propDir,
00046 beamSpot, useRefittedState, constructTsosWithErrors );
00047 }
00048 else
00049 {
00050 theValidityFlag = this->construct( trajectoryState, recHits, magField, materialEffects, propDir,
00051 beamSpot, useRefittedState, constructTsosWithErrors );
00052 }
00053 }
00054
00055
00056 TwoBodyDecayTrajectory::TwoBodyDecayTrajectory( void )
00057 : ReferenceTrajectoryBase( 0, 0, 0, 0)
00058 {}
00059
00060
00061 bool TwoBodyDecayTrajectory::construct( const TwoBodyDecayTrajectoryState& state,
00062 const ConstRecHitCollection& recHits,
00063 const MagneticField* field,
00064 MaterialEffects materialEffects,
00065 PropagationDirection propDir,
00066 const reco::BeamSpot &beamSpot,
00067 bool useRefittedState,
00068 bool constructTsosWithErrors )
00069 {
00070 const TwoBodyDecayTrajectoryState::TsosContainer& tsos = state.trajectoryStates( useRefittedState );
00071 const TwoBodyDecayTrajectoryState::Derivatives& deriv = state.derivatives();
00072 double mass = state.particleMass();
00073
00074
00075
00076
00077
00078
00079 ReferenceTrajectory trajectory1( tsos.first, recHits.first, false, field, materialEffects,
00080 propDir, mass, false, beamSpot);
00081
00082
00083 if ( !trajectory1.isValid() ) return false;
00084
00085 int nLocal = deriv.first.num_row();
00086 int nTbd = deriv.first.num_col();
00087 unsigned int nHitMeas1 = trajectory1.numberOfHitMeas();
00088 unsigned int nMsMeas1 = trajectory1.numberOfMsMeas();
00089 unsigned int nPar1 = trajectory1.numberOfPar();
00090 unsigned int nMsPar1 = trajectory1.numberOfMsPar();
00091
00092
00093 AlgebraicMatrix fullDeriv1 = trajectory1.derivatives().sub(1,nHitMeas1+nMsMeas1,1,nLocal) * trajectory1.localToTrajectory() * deriv.first;
00094
00095
00096
00097
00098
00099 ReferenceTrajectory trajectory2( tsos.second, recHits.second, false, field, materialEffects,
00100 propDir, mass, false, beamSpot );
00101
00102 if ( !trajectory2.isValid() ) return false;
00103
00104 unsigned int nHitMeas2 = trajectory2.numberOfHitMeas();
00105 unsigned int nMsMeas2 = trajectory2.numberOfMsMeas();
00106 unsigned int nPar2 = trajectory2.numberOfPar();
00107 unsigned int nMsPar2 = trajectory2.numberOfMsPar();
00108
00109 AlgebraicMatrix fullDeriv2 = trajectory2.derivatives().sub(1,nHitMeas2+nMsMeas2,1,nLocal) * trajectory2.localToTrajectory() * deriv.second;
00110
00111
00112
00113
00114
00115 theNumberOfRecHits.first = recHits.first.size();
00116 theNumberOfRecHits.second = recHits.second.size();
00117
00118 theNumberOfHits = trajectory1.numberOfHits() + trajectory2.numberOfHits();
00119 theNumberOfPars = nPar1 + nPar2;
00120 theNumberOfMsPars = nMsPar1 + nMsPar2;
00121 theNumberOfMsMeas = nMsMeas1 + nMsMeas2;
00122
00123
00124 int rowOffset = 1;
00125 int colOffset = 1;
00126 theDerivatives.sub( rowOffset, colOffset, fullDeriv1.sub( 1, nHitMeas1, 1, nTbd ) );
00127 colOffset += nTbd;
00128 theDerivatives.sub( rowOffset, colOffset, trajectory1.derivatives().sub( 1, nHitMeas1, nLocal + 1, nPar1 + nMsPar1 ) );
00129
00130 rowOffset += nHitMeas1;
00131 colOffset = 1;
00132 theDerivatives.sub( rowOffset, colOffset, fullDeriv2.sub( 1, nHitMeas2, 1, nTbd ) );
00133 colOffset += (nPar1 + nMsPar1 + nTbd - nLocal);
00134 theDerivatives.sub( rowOffset, colOffset, trajectory2.derivatives().sub( 1, nHitMeas2, nLocal + 1, nPar2 + nMsPar2 ) );
00135
00136 rowOffset += nHitMeas2;
00137 colOffset = 1;
00138 theDerivatives.sub( rowOffset, colOffset, fullDeriv1.sub(nHitMeas1 + 1, nHitMeas1 + nMsMeas1, 1, nTbd ) );
00139 colOffset += nTbd;
00140 theDerivatives.sub( rowOffset, colOffset, trajectory1.derivatives().sub(nHitMeas1 + 1, nHitMeas1 + nMsMeas1, nLocal + 1, nPar1 + nMsPar1 ) );
00141
00142 rowOffset += nMsMeas1;
00143 colOffset = 1;
00144 theDerivatives.sub( rowOffset, colOffset, fullDeriv2.sub(nHitMeas2 + 1, nHitMeas2 + nMsMeas2, 1, nTbd ) );
00145 colOffset += (nPar1 + nMsPar1 + nTbd - nLocal);
00146 theDerivatives.sub( rowOffset, colOffset, trajectory2.derivatives().sub(nHitMeas2 + 1, nHitMeas2 + nMsMeas2, nLocal + 1, nPar2 + nMsPar2 ) );
00147
00148 theMeasurements.sub( 1, trajectory1.measurements().sub( 1, nHitMeas1 ) );
00149 theMeasurements.sub( nHitMeas1 + 1, trajectory2.measurements().sub( 1, nHitMeas2 ) );
00150 theMeasurements.sub( nHitMeas1 + nHitMeas2 + 1, trajectory1.measurements().sub(nHitMeas1 + 1, nHitMeas1 + nMsMeas1 ) );
00151 theMeasurements.sub( nHitMeas1 + nHitMeas2 + nMsMeas1 + 1, trajectory2.measurements().sub(nHitMeas2 + 1, nHitMeas2 + nMsMeas2 ) );
00152
00153 theMeasurementsCov.sub( 1, trajectory1.measurementErrors().sub( 1, nHitMeas1 ) );
00154 theMeasurementsCov.sub( nHitMeas1 + 1, trajectory2.measurementErrors().sub( 1, nHitMeas2 ) );
00155 theMeasurementsCov.sub( nHitMeas1 + nHitMeas2 + 1, trajectory1.measurementErrors().sub(nHitMeas1 + 1, nHitMeas1 + nMsMeas1 ) );
00156 theMeasurementsCov.sub( nHitMeas1 + nHitMeas2 + nMsMeas1 + 1, trajectory2.measurementErrors().sub(nHitMeas2 + 1, nHitMeas2 + nMsMeas2 ) );
00157
00158 theTrajectoryPositions.sub( 1, trajectory1.trajectoryPositions() );
00159 theTrajectoryPositions.sub( nHitMeas1 + 1, trajectory2.trajectoryPositions() );
00160
00161 theTrajectoryPositionCov = state.decayParameters().covariance().similarity( theDerivatives.sub(1, nHitMeas1 + nHitMeas2, 1, 9) );
00162
00163 theParameters = state.decayParameters().parameters();
00164
00165 theRecHits.insert( theRecHits.end(), recHits.first.begin(), recHits.first.end() );
00166 theRecHits.insert( theRecHits.end(), recHits.second.begin(), recHits.second.end() );
00167
00168 if ( constructTsosWithErrors )
00169 {
00170 constructTsosVecWithErrors( trajectory1, trajectory2, field );
00171 }
00172 else
00173 {
00174 theTsosVec.insert( theTsosVec.end(),
00175 trajectory1.trajectoryStates().begin(),
00176 trajectory1.trajectoryStates().end() );
00177
00178 theTsosVec.insert( theTsosVec.end(),
00179 trajectory2.trajectoryStates().begin(),
00180 trajectory2.trajectoryStates().end() );
00181 }
00182
00183 return true;
00184 }
00185
00186
00187 void TwoBodyDecayTrajectory::constructTsosVecWithErrors( const ReferenceTrajectory& traj1,
00188 const ReferenceTrajectory& traj2,
00189 const MagneticField* field )
00190 {
00191 int iTsos = 0;
00192
00193 std::vector< TrajectoryStateOnSurface >::const_iterator itTsos;
00194
00195 for ( itTsos = traj1.trajectoryStates().begin(); itTsos != traj1.trajectoryStates().end(); itTsos++ )
00196 {
00197 constructSingleTsosWithErrors( *itTsos, iTsos, field );
00198 iTsos++;
00199 }
00200
00201 for ( itTsos = traj2.trajectoryStates().begin(); itTsos != traj2.trajectoryStates().end(); itTsos++ )
00202 {
00203 constructSingleTsosWithErrors( *itTsos, iTsos, field );
00204 iTsos++;
00205 }
00206 }
00207
00208
00209 void TwoBodyDecayTrajectory::constructSingleTsosWithErrors( const TrajectoryStateOnSurface& tsos,
00210 int iTsos,
00211 const MagneticField* field )
00212 {
00213 AlgebraicSymMatrix localErrors( 5, 0 );
00214 const double coeff = 1e-2;
00215
00216 double invP = tsos.localParameters().signedInverseMomentum();
00217 LocalVector p = tsos.localParameters().momentum();
00218
00219
00220
00221 float dpinv = coeff*( fabs(p.x()) + fabs(p.y()) + fabs(p.z()) )*invP*invP;
00222 float dxdir = coeff*( fabs(p.x()) + fabs(p.z()) )/p.z()/p.z();
00223 float dydir = coeff*( fabs(p.y()) + fabs(p.z()) )/p.z()/p.z();
00224 localErrors[0][0] = dpinv*dpinv;
00225 localErrors[1][1] = dxdir*dxdir;
00226 localErrors[2][2] = dydir*dydir;
00227
00228
00229 localErrors[3][3] = theTrajectoryPositionCov[nMeasPerHit*iTsos][nMeasPerHit*iTsos];
00230 localErrors[3][4] = theTrajectoryPositionCov[nMeasPerHit*iTsos][nMeasPerHit*iTsos+1];
00231 localErrors[4][4] = theTrajectoryPositionCov[nMeasPerHit*iTsos+1][nMeasPerHit*iTsos+1];
00232
00233
00234 theTsosVec[iTsos] = TrajectoryStateOnSurface( tsos.localParameters(),
00235 LocalTrajectoryError( localErrors ),
00236 tsos.surface(),
00237 field,
00238 tsos.surfaceSide() );
00239 }