CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoVertex/VertexTools/src/SequentialVertexFitter.cc

Go to the documentation of this file.
00001 #include "RecoVertex/VertexTools/interface/SequentialVertexFitter.h"
00002 #include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 
00005 #include <algorithm>
00006 using namespace std;
00007 using namespace reco;
00008 
00009 namespace {
00010   // FIXME
00011   // hard-coded tracker bounds
00012   // workaround while waiting for Geometry service
00013   static const float TrackerBoundsRadius = 112;
00014   static const float TrackerBoundsHalfLength = 273.5;
00015   bool insideTrackerBounds(const GlobalPoint& point) {
00016     return ((point.transverse() < TrackerBoundsRadius)
00017         && (abs(point.z()) < TrackerBoundsHalfLength));
00018   }
00019 }
00020 
00021 
00022 template <unsigned int N>
00023 SequentialVertexFitter<N>::SequentialVertexFitter(
00024   const LinearizationPointFinder & linP, 
00025   const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
00026   const AbstractLTSFactory<N> & ltsf ) : 
00027   theLinP(linP.clone()), theUpdator(updator.clone()), 
00028   theSmoother(smoother.clone()), theLTrackFactory ( ltsf.clone() )
00029 {
00030   setDefaultParameters();
00031 }
00032 
00033 template <unsigned int N>
00034 SequentialVertexFitter<N>::SequentialVertexFitter(
00035   const edm::ParameterSet& pSet, const LinearizationPointFinder & linP, 
00036   const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
00037   const AbstractLTSFactory<N> & ltsf) : 
00038   thePSet(pSet), theLinP(linP.clone()), theUpdator(updator.clone()), 
00039   theSmoother(smoother.clone()), theLTrackFactory ( ltsf.clone() )
00040 {
00041   readParameters();
00042 }
00043 
00044 
00045 template <unsigned int N>
00046 SequentialVertexFitter<N>::SequentialVertexFitter(
00047   const SequentialVertexFitter & original)
00048 {
00049   thePSet = original.parameterSet();
00050   theLinP = original.linearizationPointFinder()->clone();
00051   theUpdator = original.vertexUpdator()->clone();
00052   theSmoother = original.vertexSmoother()->clone();
00053   theMaxShift = original.maxShift();
00054   theMaxStep = original.maxStep();
00055   theLTrackFactory = original.linearizedTrackStateFactory()->clone();
00056 }
00057 
00058 
00059 template <unsigned int N>
00060 SequentialVertexFitter<N>::~SequentialVertexFitter()
00061 {
00062   delete theLinP;
00063   delete theUpdator;
00064   delete theSmoother;
00065   delete theLTrackFactory;
00066 }
00067 
00068 
00069 template <unsigned int N>
00070 void SequentialVertexFitter<N>::readParameters()
00071 {
00072   theMaxShift = thePSet.getParameter<double>("maxDistance"); //0.01
00073   theMaxStep = thePSet.getParameter<int>("maxNbrOfIterations"); //10
00074 }
00075 
00076 template <unsigned int N>
00077 void SequentialVertexFitter<N>::setDefaultParameters()
00078 {
00079   thePSet.addParameter<double>("maxDistance", 0.01);
00080   thePSet.addParameter<int>("maxNbrOfIterations", 10); //10
00081   readParameters();
00082 }
00083 
00084 template <unsigned int N>
00085 CachingVertex<N> 
00086 SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack> & tracks) const
00087 {
00088   // Linearization Point
00089   GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
00090   if (!insideTrackerBounds(linP)) linP = GlobalPoint(0,0,0);
00091 
00092   // Initial vertex state, with a very large error matrix
00093   ROOT::Math::SMatrixIdentity id;
00094   AlgebraicSymMatrix33 we(id);
00095   GlobalError error(we*10000);
00096   VertexState state(linP, error);
00097   std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
00098   return fit(vtContainer, state, false);
00099 }
00100 
00101 template <unsigned int N>
00102 CachingVertex<N> SequentialVertexFitter<N>::vertex(
00103     const std::vector<RefCountedVertexTrack> & tracks,
00104     const reco::BeamSpot & spot ) const
00105 {
00106   VertexState state(spot);
00107   return fit(tracks, state, true );
00108 }
00109 
00110 template <unsigned int N>
00111 CachingVertex<N> 
00112 SequentialVertexFitter<N>::vertex(const std::vector<RefCountedVertexTrack> & tracks) const
00113 {
00114   // Initial vertex state, with a very small weight matrix
00115   GlobalPoint linP = tracks[0]->linearizedTrack()->linearizationPoint();
00116   ROOT::Math::SMatrixIdentity id;
00117   AlgebraicSymMatrix33 we(id);
00118   GlobalError error(we*10000);
00119   VertexState state(linP, error);
00120   return fit(tracks, state, false);
00121 }
00122 
00123 
00124 // Fit vertex out of a set of RecTracks. 
00125 // Uses the specified linearization point.
00126 //
00127 template <unsigned int N>
00128 CachingVertex<N>  
00129 SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack> & tracks, 
00130                                const GlobalPoint& linPoint) const
00131 { 
00132   // Initial vertex state, with a very large error matrix
00133   ROOT::Math::SMatrixIdentity id;
00134   AlgebraicSymMatrix33 we(id);
00135   GlobalError error(we*10000);
00136   VertexState state(linPoint, error);
00137   std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
00138   return fit(vtContainer, state, false);
00139 }
00140 
00141 
00146 template <unsigned int N>
00147 CachingVertex<N> 
00148 SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack> & tracks,
00149                                const BeamSpot& beamSpot) const
00150 {
00151   VertexState beamSpotState(beamSpot);
00152   std::vector<RefCountedVertexTrack> vtContainer;
00153 
00154   if (tracks.size() > 1) {
00155     // Linearization Point search if there are more than 1 track
00156     GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
00157     if (!insideTrackerBounds(linP)) linP = GlobalPoint(0,0,0);
00158     ROOT::Math::SMatrixIdentity id;
00159     AlgebraicSymMatrix33 we(id);
00160     GlobalError error(we*10000);
00161     VertexState lpState(linP, error);
00162     vtContainer = linearizeTracks(tracks, lpState);
00163   } else {
00164     // otherwise take the beamspot position.
00165     vtContainer = linearizeTracks(tracks, beamSpotState);
00166   }
00167 
00168   return fit(vtContainer, beamSpotState, true);
00169 }
00170 
00171 
00172 // Fit vertex out of a set of RecTracks. 
00173 // Uses the position as both the linearization point AND as prior
00174 // estimate of the vertex position. The error is used for the 
00175 // weight of the prior estimate.
00176 //
00177 template <unsigned int N>
00178 CachingVertex<N> SequentialVertexFitter<N>::vertex(
00179   const std::vector<reco::TransientTrack> & tracks, 
00180   const GlobalPoint& priorPos,
00181   const GlobalError& priorError) const
00182 { 
00183   VertexState state(priorPos, priorError);
00184   std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
00185   return fit(vtContainer, state, true);
00186 }
00187 
00188 // Fit vertex out of a set of VertexTracks
00189 // Uses the position and error for the prior estimate of the vertex.
00190 // This position is not used to relinearize the tracks.
00191 //
00192 template <unsigned int N>
00193 CachingVertex<N> SequentialVertexFitter<N>::vertex(
00194   const std::vector<RefCountedVertexTrack> & tracks, 
00195   const GlobalPoint& priorPos,
00196   const GlobalError& priorError) const
00197 {
00198   VertexState state(priorPos, priorError);
00199   return fit(tracks, state, true);
00200 }
00201 
00202 
00203 // Construct a container of VertexTrack from a set of RecTracks.
00204 //
00205 template <unsigned int N>
00206 vector<typename SequentialVertexFitter<N>::RefCountedVertexTrack> 
00207 SequentialVertexFitter<N>::linearizeTracks(
00208   const std::vector<reco::TransientTrack> & tracks, 
00209   const VertexState state) const
00210 {
00211   GlobalPoint linP = state.position();
00212   std::vector<RefCountedVertexTrack> finalTracks;
00213   finalTracks.reserve(tracks.size());
00214   for(vector<reco::TransientTrack>::const_iterator i = tracks.begin(); 
00215       i != tracks.end(); i++) {
00216     RefCountedLinearizedTrackState lTrData 
00217       = theLTrackFactory->linearizedTrackState(linP, *i);
00218     RefCountedVertexTrack vTrData = theVTrackFactory.vertexTrack(lTrData,state);
00219     finalTracks.push_back(vTrData);
00220   }
00221   return finalTracks;
00222 }
00223 
00224 
00225 // Construct new a container of VertexTrack with a new linearization point
00226 // and vertex state, from an existing set of VertexTrack, from which only the 
00227 // recTracks will be used.
00228 //
00229 template <unsigned int N>
00230 vector<typename SequentialVertexFitter<N>::RefCountedVertexTrack> 
00231 SequentialVertexFitter<N>::reLinearizeTracks(
00232   const std::vector<RefCountedVertexTrack> & tracks, 
00233   const VertexState state) const
00234 {
00235 
00236   GlobalPoint linP = state.position();
00237   std::vector<RefCountedVertexTrack> finalTracks;
00238   finalTracks.reserve(tracks.size());
00239   for(typename std::vector<RefCountedVertexTrack>::const_iterator i = tracks.begin(); 
00240       i != tracks.end(); i++) {
00241     RefCountedLinearizedTrackState lTrData = 
00242           (**i).linearizedTrack()->stateWithNewLinearizationPoint(linP);
00243     //    RefCountedLinearizedTrackState lTrData = 
00244     //      theLTrackFactory->linearizedTrackState(linP, 
00245     //                              (**i).linearizedTrack()->track());
00246     RefCountedVertexTrack vTrData =
00247       theVTrackFactory.vertexTrack(lTrData,state, (**i).weight() );
00248     finalTracks.push_back(vTrData);
00249   }
00250   return finalTracks;
00251 }
00252 
00253 
00254 // The method where the vertex fit is actually done!
00255 //
00256 template <unsigned int N>
00257 CachingVertex<N> 
00258 SequentialVertexFitter<N>::fit(const std::vector<RefCountedVertexTrack> & tracks,
00259                             const VertexState priorVertex, bool withPrior ) const
00260 {
00261   std::vector<RefCountedVertexTrack> initialTracks;
00262   GlobalPoint priorVertexPosition = priorVertex.position();
00263   GlobalError priorVertexError = priorVertex.error();
00264   
00265   CachingVertex<N> returnVertex(priorVertexPosition,priorVertexError,initialTracks,0);
00266   if (withPrior) {
00267     returnVertex = CachingVertex<N>(priorVertexPosition,priorVertexError,
00268                 priorVertexPosition,priorVertexError,initialTracks,0);
00269   }
00270   CachingVertex<N> initialVertex = returnVertex;
00271   std::vector<RefCountedVertexTrack> globalVTracks = tracks;
00272 
00273   // main loop through all the VTracks
00274   bool validVertex = true;
00275   int step = 0;
00276   GlobalPoint newPosition = priorVertexPosition;
00277   GlobalPoint previousPosition;
00278   do {
00279     CachingVertex<N> fVertex = initialVertex;
00280     // make new linearized and vertex tracks for the next iteration
00281     if(step != 0) globalVTracks = reLinearizeTracks(tracks, 
00282                                         returnVertex.vertexState());
00283 
00284     // update sequentially the vertex estimate
00285     for (typename std::vector<RefCountedVertexTrack>::const_iterator i 
00286            = globalVTracks.begin(); i != globalVTracks.end(); i++) {
00287       fVertex = theUpdator->add(fVertex,*i);
00288       if (!fVertex.isValid()) break;
00289     }
00290 
00291     validVertex = fVertex.isValid();
00292     // check tracker bounds and NaN in position
00293     if (validVertex && hasNan(fVertex.position())) {
00294       LogDebug("RecoVertex/SequentialVertexFitter") 
00295          << "Fitted position is NaN.\n";
00296       validVertex = false;
00297     }
00298 
00299     if (validVertex && !insideTrackerBounds(fVertex.position())) {
00300       LogDebug("RecoVertex/SequentialVertexFitter") 
00301          << "Fitted position is out of tracker bounds.\n";
00302       validVertex = false;
00303     }
00304 
00305     if (!validVertex) {
00306       // reset initial vertex position to (0,0,0) and force new iteration 
00307       // if number of steps not exceeded
00308       ROOT::Math::SMatrixIdentity id;
00309       AlgebraicSymMatrix33 we(id);
00310       GlobalError error(we*10000);
00311       fVertex = CachingVertex<N>(GlobalPoint(0,0,0), error,
00312                               initialTracks, 0);
00313     }
00314 
00315     previousPosition = newPosition;
00316     newPosition = fVertex.position();
00317 
00318     returnVertex = fVertex;
00319     globalVTracks.clear();
00320     step++;
00321   } while ( (step != theMaxStep) &&
00322             (((previousPosition - newPosition).transverse() > theMaxShift) ||
00323                 (!validVertex) ) );
00324 
00325   if (!validVertex) {
00326     LogDebug("RecoVertex/SequentialVertexFitter") 
00327        << "Fitted position is invalid (out of tracker bounds or has NaN). Returned vertex is invalid\n";
00328     return CachingVertex<N>(); // return invalid vertex
00329   }
00330 
00331   if (step >= theMaxStep) {
00332     LogDebug("RecoVertex/SequentialVertexFitter") 
00333        << "The maximum number of steps has been exceeded. Returned vertex is invalid\n";
00334     return CachingVertex<N>(); // return invalid vertex
00335   }
00336 
00337   // smoothing
00338   returnVertex = theSmoother->smooth(returnVertex);
00339 
00340   return returnVertex;
00341 }
00342 
00343 template class SequentialVertexFitter<5>;
00344 template class SequentialVertexFitter<6>;