CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/Alignment/ReferenceTrajectories/plugins/DualTrajectoryFactory.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/ESHandle.h"
00002 #include "FWCore/Framework/interface/EventSetup.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h" 
00006 #include "MagneticField/Engine/interface/MagneticField.h"
00007 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryPlugin.h"
00008 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00009 
00010 #include <algorithm>
00011 
00012 #include "Alignment/ReferenceTrajectories/interface/DualReferenceTrajectory.h" 
00013 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryBase.h"
00014 
00016 
00017 
00018 class DualTrajectoryFactory : public TrajectoryFactoryBase
00019 {
00020 public:
00021   DualTrajectoryFactory(const edm::ParameterSet &config);
00022   virtual ~DualTrajectoryFactory();
00023 
00025   virtual const ReferenceTrajectoryCollection trajectories(const edm::EventSetup  &setup,
00026                                                            const ConstTrajTrackPairCollection &tracks,
00027                                                            const reco::BeamSpot &beamSpot) const;
00028 
00029   virtual const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup,
00030                                                            const ConstTrajTrackPairCollection &tracks,
00031                                                            const ExternalPredictionCollection &external,
00032                                                            const reco::BeamSpot &beamSpot) const;
00033 
00034   virtual DualTrajectoryFactory* clone() const { return new DualTrajectoryFactory(*this); }
00035 
00036 protected:
00037   struct DualTrajectoryInput
00038   {
00039     TrajectoryStateOnSurface refTsos;
00040     TransientTrackingRecHit::ConstRecHitContainer fwdRecHits;
00041     TransientTrackingRecHit::ConstRecHitContainer bwdRecHits;
00042   };
00043 
00044   const DualTrajectoryInput referenceStateAndRecHits(const ConstTrajTrackPair &track) const;
00045 
00046   const TrajectoryStateOnSurface propagateExternal(const TrajectoryStateOnSurface &external,
00047                                                    const Surface &surface,
00048                                                    const MagneticField *magField) const;
00049 
00050   double theMass;
00051 };
00052 
00056 
00057 DualTrajectoryFactory::DualTrajectoryFactory( const edm::ParameterSet & config ) :
00058   TrajectoryFactoryBase( config ),
00059   theMass(config.getParameter<double>("ParticleMass"))
00060 {
00061   edm::LogInfo("Alignment") << "@SUB=DualTrajectoryFactory"
00062                             << "mass: " << theMass;
00063 }
00064 
00065  
00066 DualTrajectoryFactory::~DualTrajectoryFactory( void ) {}
00067 
00068 
00069 const DualTrajectoryFactory::ReferenceTrajectoryCollection
00070 DualTrajectoryFactory::trajectories(const edm::EventSetup &setup,
00071                                     const ConstTrajTrackPairCollection &tracks,
00072                                     const reco::BeamSpot &beamSpot) const
00073 {
00074   ReferenceTrajectoryCollection trajectories;
00075 
00076   edm::ESHandle< MagneticField > magneticField;
00077   setup.get< IdealMagneticFieldRecord >().get( magneticField );
00078   if (magneticField->inTesla(GlobalPoint(0.,0.,0.)).mag2() < 1.e-6) {
00079     edm::LogWarning("Alignment") << "@SUB=DualTrajectoryFactory::trajectories"
00080                                  << "B-field in z is " << magneticField->inTesla(GlobalPoint(0.,0.,0.)).z()
00081                                  << ": You should probably use the DualBzeroTrajectoryFactory\n" 
00082                                  << "or fix this code to switch automatically as the ReferenceTrajectoryFactory does.";
00083   }
00084 
00085   ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
00086 
00087   while ( itTracks != tracks.end() )
00088   { 
00089     const DualTrajectoryInput input = this->referenceStateAndRecHits( *itTracks );
00090     // Check input: If all hits were rejected, the TSOS is initialized as invalid.
00091     if ( input.refTsos.isValid() )
00092     {
00093       ReferenceTrajectoryPtr ptr( new DualReferenceTrajectory( input.refTsos,
00094                                                                input.fwdRecHits,
00095                                                                input.bwdRecHits,
00096                                                                magneticField.product(),
00097                                                                materialEffects(),
00098                                                                propagationDirection(),
00099                                                                theMass,
00100                                                                theUseBeamSpot, beamSpot ) );
00101       trajectories.push_back( ptr );
00102     }
00103 
00104     ++itTracks;
00105   }
00106 
00107   return trajectories;
00108 }
00109 
00110 const DualTrajectoryFactory::ReferenceTrajectoryCollection
00111 DualTrajectoryFactory::trajectories(const edm::EventSetup &setup,
00112                                     const ConstTrajTrackPairCollection &tracks,
00113                                     const ExternalPredictionCollection &external,
00114                                     const reco::BeamSpot &beamSpot) const
00115 {
00116   ReferenceTrajectoryCollection trajectories;
00117 
00118   if ( tracks.size() != external.size() )
00119   {
00120     edm::LogInfo("ReferenceTrajectories") << "@SUB=DualTrajectoryFactory::trajectories"
00121                                           << "Inconsistent input:\n"
00122                                           << "\tnumber of tracks = " << tracks.size()
00123                                           << "\tnumber of external predictions = " << external.size();
00124     return trajectories;
00125   }
00126 
00127   edm::ESHandle< MagneticField > magneticField;
00128   setup.get< IdealMagneticFieldRecord >().get( magneticField );
00129   if (magneticField->inTesla(GlobalPoint(0.,0.,0.)).mag2() < 1.e-6) {
00130     edm::LogWarning("Alignment") << "@SUB=DualTrajectoryFactory::trajectories"
00131                                  << "B-field in z is " << magneticField->inTesla(GlobalPoint(0.,0.,0.)).z()
00132                                  << ": You should probably use the DualBzeroTrajectoryFactory\n" 
00133                                  << "or fix this code to switch automatically as the ReferenceTrajectoryFactory does.";
00134   }
00135 
00136   ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
00137   ExternalPredictionCollection::const_iterator itExternal = external.begin();
00138 
00139   while ( itTracks != tracks.end() )
00140   {
00141     const DualTrajectoryInput input = referenceStateAndRecHits( *itTracks );
00142     // Check input: If all hits were rejected, the TSOS is initialized as invalid.
00143     if ( input.refTsos.isValid() )
00144     {
00145       if ( (*itExternal).isValid() )
00146       {
00147         TrajectoryStateOnSurface propExternal =
00148           propagateExternal( *itExternal, input.refTsos.surface(), magneticField.product() );
00149 
00150         if ( !propExternal.isValid() ) continue;
00151 
00152         // set the flag for reversing the RecHits to false, since they are already in the correct order.
00153         ReferenceTrajectoryPtr ptr( new DualReferenceTrajectory( propExternal,
00154                                                                  input.fwdRecHits,
00155                                                                  input.bwdRecHits,
00156                                                                  magneticField.product(),
00157                                                                  materialEffects(),
00158                                                                  propagationDirection(),
00159                                                                  theMass, 
00160                                                                  theUseBeamSpot, beamSpot ) );
00161 
00162         AlgebraicSymMatrix externalParamErrors( asHepMatrix<5>( propExternal.localError().matrix() ) );
00163         ptr->setParameterErrors( externalParamErrors );
00164         trajectories.push_back( ptr );
00165       }
00166       else
00167       {
00168         ReferenceTrajectoryPtr ptr( new DualReferenceTrajectory( input.refTsos,
00169                                                                  input.fwdRecHits,
00170                                                                  input.bwdRecHits,
00171                                                                  magneticField.product(),
00172                                                                  materialEffects(),
00173                                                                  propagationDirection(),
00174                                                                  theMass, 
00175                                                                  theUseBeamSpot, beamSpot ) );
00176         trajectories.push_back( ptr );
00177       }
00178     }
00179 
00180     ++itTracks;
00181     ++itExternal;
00182   }
00183 
00184   return trajectories;
00185 }
00186 
00187 
00188 const DualTrajectoryFactory::DualTrajectoryInput
00189 DualTrajectoryFactory::referenceStateAndRecHits( const ConstTrajTrackPair& track ) const
00190 {
00191   DualTrajectoryInput input;
00192  
00193   // get the trajectory measurements in the correct order, i.e. reverse if needed
00194   Trajectory::DataContainer allTrajMeas = this->orderedTrajectoryMeasurements( *track.first );
00195   Trajectory::DataContainer usedTrajMeas;
00196   Trajectory::DataContainer::iterator itM;
00197   // get all relevant trajectory measurements
00198   for ( itM = allTrajMeas.begin(); itM != allTrajMeas.end(); itM++ )
00199   {
00200     if ( useRecHit( ( *itM ).recHit() ) ) usedTrajMeas.push_back( *itM );
00201   }
00202 
00203   unsigned int iMeas = 0;
00204   unsigned int nMeas = usedTrajMeas.size();
00205   unsigned int nRefStateMeas = nMeas/2;
00206   // get the valid RecHits
00207   for ( itM = usedTrajMeas.begin(); itM != usedTrajMeas.end(); itM++, iMeas++ )
00208   {
00209     TransientTrackingRecHit::ConstRecHitPointer aRecHit = ( *itM ).recHit();
00210 
00211     if ( iMeas < nRefStateMeas ) {
00212       input.bwdRecHits.push_back( aRecHit );
00213     } else if ( iMeas > nRefStateMeas ) {
00214       input.fwdRecHits.push_back( aRecHit );
00215     } else { // iMeas == nRefStateMeas
00216       if ( ( *itM ).updatedState().isValid() )
00217       {
00218         input.refTsos = ( *itM ).updatedState();
00219         input.bwdRecHits.push_back( aRecHit );
00220         input.fwdRecHits.push_back( aRecHit );
00221       } else {
00222         // if the tsos of the middle hit is not valid, try the next one ...
00223         nRefStateMeas++;
00224         input.bwdRecHits.push_back( aRecHit );
00225       }
00226     }
00227   }
00228 
00229   // bring input.fwdRecHits into correct order
00230   std::reverse( input.bwdRecHits.begin(), input.bwdRecHits.end() );
00231 
00232   return input;
00233 }
00234 
00235 const TrajectoryStateOnSurface
00236 DualTrajectoryFactory::propagateExternal( const TrajectoryStateOnSurface& external,
00237                                           const Surface& surface,
00238                                           const MagneticField* magField ) const
00239 {
00240   AnalyticalPropagator propagator( magField, anyDirection );
00241   const std::pair< TrajectoryStateOnSurface, double > tsosWithPath =
00242     propagator.propagateWithPath( external, surface );
00243   return tsosWithPath.first;
00244 }
00245 
00246 
00247 DEFINE_EDM_PLUGIN( TrajectoryFactoryPlugin, DualTrajectoryFactory, "DualTrajectoryFactory" );