CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/Alignment/ReferenceTrajectories/src/DualReferenceTrajectory.cc

Go to the documentation of this file.
00001 
00002 #include "Alignment/ReferenceTrajectories/interface/DualReferenceTrajectory.h"
00003 
00004 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00005 
00006 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h" 
00007 #include "DataFormats/TrajectoryState/interface/LocalTrajectoryParameters.h"
00008 
00009 #include "Alignment/ReferenceTrajectories/interface/ReferenceTrajectory.h"
00010 
00011 DualReferenceTrajectory::DualReferenceTrajectory( const TrajectoryStateOnSurface &referenceTsos,
00012                                                   const ConstRecHitContainer &forwardRecHits,
00013                                                   const ConstRecHitContainer &backwardRecHits,
00014                                                   const MagneticField *magField,
00015                                                   MaterialEffects materialEffects,
00016                                                   PropagationDirection propDir,
00017                                                   double mass, 
00018                                                   bool useBeamSpot,
00019                                                   const reco::BeamSpot &beamSpot)
00020   : ReferenceTrajectoryBase(referenceTsos.localParameters().mixedFormatVector().kSize,
00021                             numberOfUsedRecHits(forwardRecHits) + numberOfUsedRecHits(backwardRecHits) - 1,
00022                             0, 0)
00023 {
00024   theValidityFlag = this->construct( referenceTsos, forwardRecHits, backwardRecHits,
00025                                      mass, materialEffects, propDir, magField, useBeamSpot, beamSpot );
00026 }
00027 
00028 
00029 
00030 DualReferenceTrajectory::DualReferenceTrajectory(unsigned int nPar, unsigned int nHits)
00031   : ReferenceTrajectoryBase(nPar, nHits, 0, 0)
00032 {}
00033 
00034 
00035 bool DualReferenceTrajectory::construct( const TrajectoryStateOnSurface &refTsos, 
00036                                          const ConstRecHitContainer &forwardRecHits,
00037                                          const ConstRecHitContainer &backwardRecHits,
00038                                          double mass, MaterialEffects materialEffects,
00039                                          const PropagationDirection propDir,
00040                                          const MagneticField *magField,
00041                                          bool useBeamSpot,
00042                                          const reco::BeamSpot &beamSpot)
00043 {
00044   if (materialEffects >= breakPoints)  throw cms::Exception("BadConfig")
00045     << "[DualReferenceTrajectory::construct] Wrong MaterialEffects: " << materialEffects;
00046     
00047   ReferenceTrajectoryBase* fwdTraj = construct(refTsos, forwardRecHits,
00048                                                mass, materialEffects,
00049                                                propDir, magField, useBeamSpot, beamSpot);
00050 
00051   ReferenceTrajectoryBase* bwdTraj = construct(refTsos, backwardRecHits,
00052                                                mass, materialEffects,
00053                                                oppositeDirection(propDir), magField, false, beamSpot);
00054 
00055   if ( !( fwdTraj->isValid() && bwdTraj->isValid() ) )
00056   {
00057     delete fwdTraj;
00058     delete bwdTraj;
00059     return false;
00060   }
00061 
00062   //
00063   // Combine both reference trajactories to a dual reference trajectory
00064   //
00065 
00066   const std::vector<TrajectoryStateOnSurface>& fwdTsosVec = fwdTraj->trajectoryStates();
00067   const std::vector<TrajectoryStateOnSurface>& bwdTsosVec = bwdTraj->trajectoryStates();
00068   theTsosVec.insert( theTsosVec.end(), fwdTsosVec.begin(), fwdTsosVec.end() );
00069   theTsosVec.insert( theTsosVec.end(), ++bwdTsosVec.begin(), bwdTsosVec.end() );
00070 
00071   const ConstRecHitContainer &fwdRecHits = fwdTraj->recHits();
00072   const ConstRecHitContainer &bwdRecHits = bwdTraj->recHits();
00073   theRecHits.insert( theRecHits.end(), fwdRecHits.begin(), fwdRecHits.end() );
00074   theRecHits.insert( theRecHits.end(), ++bwdRecHits.begin(), bwdRecHits.end() );
00075 
00076   theParameters = extractParameters( refTsos );
00077   
00078   unsigned int nParam   = theNumberOfPars;
00079   unsigned int nFwdMeas = fwdTraj->numberOfHitMeas();
00080   unsigned int nBwdMeas = bwdTraj->numberOfHitMeas();
00081   unsigned int nFwdBP   = fwdTraj->numberOfVirtualMeas();
00082   unsigned int nBwdBP   = bwdTraj->numberOfVirtualMeas();
00083   unsigned int nMeas    = nFwdMeas+nBwdMeas-nMeasPerHit; 
00084        
00085   theMeasurements.sub( 1, fwdTraj->measurements().sub( 1, nFwdMeas ) );
00086   theMeasurements.sub( nFwdMeas+1, bwdTraj->measurements().sub( nMeasPerHit+1, nBwdMeas ) );
00087     
00088   theMeasurementsCov.sub( 1, fwdTraj->measurementErrors().sub( 1, nFwdMeas ) );
00089   theMeasurementsCov.sub( nFwdMeas+1, bwdTraj->measurementErrors().sub( nMeasPerHit+1, nBwdMeas ) );
00090   
00091   theTrajectoryPositions.sub( 1, fwdTraj->trajectoryPositions() );
00092   theTrajectoryPositions.sub( nFwdMeas+1, bwdTraj->trajectoryPositions().sub( nMeasPerHit+1, nBwdMeas ) );
00093 
00094   theTrajectoryPositionCov.sub( 1, fwdTraj->trajectoryPositionErrors() );
00095   theTrajectoryPositionCov.sub( nFwdMeas+1, bwdTraj->trajectoryPositionErrors().sub( nMeasPerHit+1, nBwdMeas ) );
00096 
00097   theDerivatives.sub( 1, 1, fwdTraj->derivatives().sub( 1, nFwdMeas, 1, nParam ) );
00098   theDerivatives.sub( nFwdMeas+1, 1, bwdTraj->derivatives().sub( nMeasPerHit+1, nBwdMeas, 1, nParam ) );
00099   
00100 // for the break points 
00101 // DUAL with break points makes no sense: (MS) correlations between the two parts are lost ! 
00102   if (nFwdBP>0 )
00103   {
00104     theMeasurements.sub( nMeas+1, fwdTraj->measurements().sub( nFwdMeas+1, nFwdMeas+nFwdBP ) );
00105     theMeasurementsCov.sub( nMeas+1, fwdTraj->measurementErrors().sub( nFwdMeas+1, nFwdMeas+nFwdBP ) );  
00106     theDerivatives.sub( 1, nParam+1, fwdTraj->derivatives().sub( 1, nFwdMeas, nParam+1, nParam+nFwdBP ) );
00107     theDerivatives.sub( nMeas+1, nParam+1, fwdTraj->derivatives().sub( nFwdMeas+1, nFwdMeas+nFwdBP, nParam+1, nParam+nFwdBP ) );
00108   }  
00109   if (nBwdBP>0 )
00110   {
00111     theMeasurements.sub( nMeas+nFwdBP+1, bwdTraj->measurements().sub( nBwdMeas+1, nBwdMeas+nBwdBP ) );  
00112     theMeasurementsCov.sub( nMeas+nFwdBP+1, bwdTraj->measurementErrors().sub( nBwdMeas+1, nBwdMeas+nBwdBP ) );  
00113     theDerivatives.sub( nFwdMeas+1, nParam+nFwdBP+1, bwdTraj->derivatives().sub( nMeasPerHit+1, nBwdMeas, nParam+1, nParam+nBwdBP ) );
00114     theDerivatives.sub( nMeas+nFwdBP+1, nParam+nFwdBP+1, bwdTraj->derivatives().sub( nBwdMeas+1, nBwdMeas+nBwdBP, nParam+1, nParam+nBwdBP ) );
00115   }
00116     
00117   delete fwdTraj;
00118   delete bwdTraj;
00119 
00120   return true; 
00121 }
00122 
00123 
00124 ReferenceTrajectory*
00125 DualReferenceTrajectory::construct(const TrajectoryStateOnSurface &referenceTsos, 
00126                                    const ConstRecHitContainer &recHits,
00127                                    double mass, MaterialEffects materialEffects,
00128                                    const PropagationDirection propDir,
00129                                    const MagneticField *magField,
00130                                    bool useBeamSpot,
00131                                    const reco::BeamSpot &beamSpot) const
00132 {
00133   if (materialEffects >= breakPoints)  throw cms::Exception("BadConfig")
00134     << "[DualReferenceTrajectory::construct] Wrong MaterialEffects: " << materialEffects;
00135   
00136   return new ReferenceTrajectory(referenceTsos, recHits, false,
00137                                  magField, materialEffects, propDir, mass,
00138                                  useBeamSpot, beamSpot);
00139 }
00140 
00141 
00142 AlgebraicVector
00143 DualReferenceTrajectory::extractParameters(const TrajectoryStateOnSurface &referenceTsos) const
00144 {
00145   return asHepVector<5>( referenceTsos.localParameters().mixedFormatVector() );
00146 }