CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/Alignment/ReferenceTrajectories/plugins/DualKalmanFactory.cc

Go to the documentation of this file.
00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 
00010 
00011 
00012 
00013 
00014 
00015 
00016 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryBase.h"
00017 
00018 #include "FWCore/Framework/interface/ESHandle.h"
00019 #include "FWCore/Framework/interface/EventSetup.h"
00020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h" 
00023 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryPlugin.h"
00024 // #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00025 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00026 
00027 #include "Alignment/ReferenceTrajectories/interface/DualKalmanTrajectory.h" 
00028 
00029 #include <algorithm>
00030 
00031 
00032 class DualKalmanFactory : public TrajectoryFactoryBase
00033 {
00034 
00035 public:
00036 
00037   DualKalmanFactory( const edm::ParameterSet & config );
00038   virtual ~DualKalmanFactory();
00039 
00041   virtual const ReferenceTrajectoryCollection 
00042   trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks,
00043                const reco::BeamSpot &beamSpot) const;
00044 
00045   virtual const ReferenceTrajectoryCollection 
00046   trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks,
00047                const ExternalPredictionCollection &external, const reco::BeamSpot &beamSpot) const;
00048 
00049   virtual DualKalmanFactory* clone() const { return new DualKalmanFactory(*this); }
00050 
00051 protected:
00052 
00053   struct DualKalmanInput 
00054   {
00055     TrajectoryStateOnSurface refTsos;
00056     Trajectory::DataContainer trajMeasurements;
00057     std::vector<unsigned int> fwdRecHitNums;
00058     std::vector<unsigned int> bwdRecHitNums;
00059   };
00060 
00061   const DualKalmanInput referenceStateAndRecHits(const ConstTrajTrackPair &track) const;
00062 
00063 //   const TrajectoryStateOnSurface propagateExternal(const TrajectoryStateOnSurface &external,
00064 //                                                 const Surface &surface,
00065 //                                                 const MagneticField *magField) const;
00066   
00067   const double theMass;
00068   const int theResidMethod;
00069 };
00070 
00071 
00072 //-----------------------------------------------------------------------------------------------
00073 //-----------------------------------------------------------------------------------------------
00074 //-----------------------------------------------------------------------------------------------
00075 DualKalmanFactory::DualKalmanFactory(const edm::ParameterSet &config) 
00076   : TrajectoryFactoryBase(config), theMass(config.getParameter<double>("ParticleMass")),
00077     theResidMethod(config.getParameter<int>("ResidualMethod"))
00078 {
00079   // Since theResidMethod is passed to DualKalmanTrajectory, valid values are checked there.
00080   edm::LogInfo("Alignment") << "@SUB=DualKalmanFactory" << "Factory created.";
00081 }
00082 
00083  
00084 //-----------------------------------------------------------------------------------------------
00085 DualKalmanFactory::~DualKalmanFactory() {}
00086 
00087 
00088 //-----------------------------------------------------------------------------------------------
00089 const DualKalmanFactory::ReferenceTrajectoryCollection
00090 DualKalmanFactory::trajectories(const edm::EventSetup &setup,
00091                                 const ConstTrajTrackPairCollection &tracks,
00092                                 const reco::BeamSpot &beamSpot) const
00093 {
00094   ReferenceTrajectoryCollection trajectories;
00095 
00096   edm::ESHandle<MagneticField> magneticField;
00097   setup.get<IdealMagneticFieldRecord>().get(magneticField);
00098 
00099   ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
00100 
00101   while (itTracks != tracks.end()) { 
00102     const DualKalmanInput input = this->referenceStateAndRecHits(*itTracks);
00103     // Check input: If all hits were rejected, the TSOS is initialized as invalid.
00104     if (input.refTsos.isValid()) {
00105       ReferenceTrajectoryPtr ptr(new DualKalmanTrajectory(input.trajMeasurements,
00106                                                           input.refTsos,
00107                                                           input.fwdRecHitNums,
00108                                                           input.bwdRecHitNums,
00109                                                           magneticField.product(),
00110                                                           this->materialEffects(),
00111                                                           this->propagationDirection(),
00112                                                           theMass, theUseBeamSpot, beamSpot,
00113                                                           theResidMethod));
00114       trajectories.push_back(ptr);
00115     }
00116     ++itTracks;
00117   }
00118 
00119   return trajectories;
00120 }
00121 
00122 //-----------------------------------------------------------------------------------------------
00123 const DualKalmanFactory::ReferenceTrajectoryCollection
00124 DualKalmanFactory::trajectories(const edm::EventSetup &setup,
00125                                 const ConstTrajTrackPairCollection &tracks,
00126                                 const ExternalPredictionCollection &external,
00127                                 const reco::BeamSpot &beamSpot) const
00128 {
00129   ReferenceTrajectoryCollection trajectories;
00130 
00131   edm::LogError("Alignment") << "@SUB=DualKalmanFactory::trajectories" 
00132                              << "Not implemented with ExternalPrediction.";
00133   return trajectories;
00134 }
00135 
00136 
00137 //-----------------------------------------------------------------------------------------------
00138 const DualKalmanFactory::DualKalmanInput
00139 DualKalmanFactory::referenceStateAndRecHits(const ConstTrajTrackPair &track) const
00140 {
00141   // Same idea as in DualTrajectoryFactory::referenceStateAndRecHits(..):
00142   // Split trajectory in the middle, take middle as reference and provide first
00143   // and second half of hits, each starting from this middle hit.
00144   // In contrast to DualTrajectoryFactory we deal here with indices and not rechits directly
00145   // to be able to get measurements and uncertainties later from the Trajectory that is
00146   // provided by the Kalman track fit.
00147 
00148   DualKalmanInput input;
00149  
00150   // get the trajectory measurements in the correct order, i.e. reverse if needed
00151   input.trajMeasurements = this->orderedTrajectoryMeasurements(*track.first);
00152 
00153   // get indices of relevant trajectory measurements to find middle of them
00154   std::vector<unsigned int> usedTrajMeasNums;
00155   for (unsigned int iM = 0; iM < input.trajMeasurements.size(); ++iM) {
00156     if (this->useRecHit(input.trajMeasurements[iM].recHit())) usedTrajMeasNums.push_back(iM);
00157   }
00158   unsigned int nRefStateMeas = usedTrajMeasNums.size()/2;
00159 
00160   // get the valid RecHits numbers
00161   for (unsigned int iMeas = 0; iMeas < usedTrajMeasNums.size(); ++iMeas) {
00162     if (iMeas < nRefStateMeas) {
00163       input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00164     } else if (iMeas > nRefStateMeas) {
00165       input.fwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00166     } else { // iMeas == nRefStateMeas
00167       if (input.trajMeasurements[usedTrajMeasNums[iMeas]].updatedState().isValid()) {
00168         input.refTsos = input.trajMeasurements[usedTrajMeasNums[iMeas]].updatedState();
00169         input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00170         input.fwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00171       } else {
00172         // if the hit/tsos of the middle hit is not valid, try the next one...
00173         ++nRefStateMeas; // but keep hit if only TSOS bad
00174         input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00175       }
00176     }
00177   }
00178 
00179   // bring input.fwdRecHits into correct order
00180   std::reverse(input.bwdRecHitNums.begin(), input.bwdRecHitNums.end());
00181 
00182   return input;
00183 }
00184 
00185 // //-----------------------------------------------------------------------------------------------
00186 // const TrajectoryStateOnSurface
00187 // DualKalmanFactory::propagateExternal(const TrajectoryStateOnSurface &external,
00188 //                                   const Surface &surface,
00189 //                                   const MagneticField *magField) const
00190 // {
00191 //   AnalyticalPropagator propagator(magField, anyDirection);
00192 //   const std::pair<TrajectoryStateOnSurface, double> tsosWithPath =
00193 //     propagator.propagateWithPath(external, surface);
00194 //   return tsosWithPath.first;
00195 //}
00196 
00197 
00198 DEFINE_EDM_PLUGIN( TrajectoryFactoryPlugin, DualKalmanFactory, "DualKalmanFactory" );