CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/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   // copy paste from DualTrajectoryFactory:
00136   if (tracks.size() != external.size()) {
00137     edm::LogInfo("ReferenceTrajectories") << "@SUB=DualKalmanFactory::trajectories"
00138                                           << "Inconsistent input:\n"
00139                                           << "\tnumber of tracks = " << tracks.size()
00140                                           << "\tnumber of external predictions = " << external.size();
00141     return trajectories;
00142   }
00143   /*
00144   edm::ESHandle<MagneticField> magneticField;
00145   setup.get<IdealMagneticFieldRecord>().get(magneticField);
00146 
00147   ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
00148   ExternalPredictionCollection::const_iterator itExternal = external.begin();
00149 
00150   while (itTracks != tracks.end()) {
00151     const DualKalmanInput input = this->referenceStateAndRecHits(*itTracks);
00152     // Check input: If all hits were rejected, the TSOS is initialized as invalid.
00153     if (input.refTsos.isValid()) {
00154       if ((*itExternal).isValid()) {
00155         TrajectoryStateOnSurface propExternal =
00156           this->propagateExternal(*itExternal, input.refTsos.surface(), magneticField.product());
00157         if (!propExternal.isValid()) continue;
00158         
00159         // set the flag for reversing the RecHits to false, since they are already in the correct order.
00160         ReferenceTrajectoryPtr ptr(new DualKalmanTrajectory(propExternal,
00161                                                             input.fwdRecHits,
00162                                                             input.bwdRecHits,
00163                                                             magneticField.product(),
00164                                                             materialEffects(),
00165                                                             propagationDirection(),
00166                                                             theMass, theResidMethod));
00167         
00168         AlgebraicSymMatrix externalParamErrors(asHepMatrix<5>(propExternal.localError().matrix()));
00169         ptr->setParameterErrors(externalParamErrors);
00170         trajectories.push_back(ptr);
00171       } else {
00172         ReferenceTrajectoryPtr ptr(new DualKalmanTrajectory(input.refTsos,
00173                                                             input.fwdRecHits,
00174                                                             input.bwdRecHits,
00175                                                             magneticField.product(),
00176                                                             materialEffects(),
00177                                                             propagationDirection(),
00178                                                             theMass, theResidMethod));
00179         trajectories.push_back(ptr);
00180       }
00181     }
00182 
00183     ++itTracks;
00184     ++itExternal;
00185   }
00186   */
00187 
00188   return trajectories;
00189 }
00190 
00191 
00192 //-----------------------------------------------------------------------------------------------
00193 const DualKalmanFactory::DualKalmanInput
00194 DualKalmanFactory::referenceStateAndRecHits(const ConstTrajTrackPair &track) const
00195 {
00196   // Same idea as in DualTrajectoryFactory::referenceStateAndRecHits(..):
00197   // Split trajectory in the middle, take middle as reference and provide first
00198   // and second half of hits, each starting from this middle hit.
00199   // In contrast to DualTrajectoryFactory we deal here with indices and not rechits directly
00200   // to be able to get measurements and uncertainties later from the Trajectory that is
00201   // provided by the Kalman track fit.
00202 
00203   DualKalmanInput input;
00204  
00205   // get the trajectory measurements in the correct order, i.e. reverse if needed
00206   input.trajMeasurements = this->orderedTrajectoryMeasurements(*track.first);
00207 
00208   // get indices of relevant trajectory measurements to find middle of them
00209   std::vector<unsigned int> usedTrajMeasNums;
00210   for (unsigned int iM = 0; iM < input.trajMeasurements.size(); ++iM) {
00211     if (this->useRecHit(input.trajMeasurements[iM].recHit())) usedTrajMeasNums.push_back(iM);
00212   }
00213   unsigned int nRefStateMeas = usedTrajMeasNums.size()/2;
00214 
00215   // get the valid RecHits numbers
00216   for (unsigned int iMeas = 0; iMeas < usedTrajMeasNums.size(); ++iMeas) {
00217     if (iMeas < nRefStateMeas) {
00218       input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00219     } else if (iMeas > nRefStateMeas) {
00220       input.fwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00221     } else { // iMeas == nRefStateMeas
00222       if (input.trajMeasurements[usedTrajMeasNums[iMeas]].updatedState().isValid()) {
00223         input.refTsos = input.trajMeasurements[usedTrajMeasNums[iMeas]].updatedState();
00224         input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00225         input.fwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00226       } else {
00227         // if the hit/tsos of the middle hit is not valid, try the next one...
00228         ++nRefStateMeas; // but keep hit if only TSOS bad
00229         input.bwdRecHitNums.push_back(usedTrajMeasNums[iMeas]);
00230       }
00231     }
00232   }
00233 
00234   // bring input.fwdRecHits into correct order
00235   std::reverse(input.bwdRecHitNums.begin(), input.bwdRecHitNums.end());
00236 
00237   return input;
00238 }
00239 
00240 // //-----------------------------------------------------------------------------------------------
00241 // const TrajectoryStateOnSurface
00242 // DualKalmanFactory::propagateExternal(const TrajectoryStateOnSurface &external,
00243 //                                   const Surface &surface,
00244 //                                   const MagneticField *magField) const
00245 // {
00246 //   AnalyticalPropagator propagator(magField, anyDirection);
00247 //   const std::pair<TrajectoryStateOnSurface, double> tsosWithPath =
00248 //     propagator.propagateWithPath(external, surface);
00249 //   return tsosWithPath.first;
00250 //}
00251 
00252 
00253 DEFINE_EDM_PLUGIN( TrajectoryFactoryPlugin, DualKalmanFactory, "DualKalmanFactory" );