CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/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   AlgebraicSymMatrix we(3,1);
00094   GlobalError error(we*10000);
00095   VertexState state(linP, error);
00096   std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
00097   return fit(vtContainer, state, false);
00098 }
00099 
00100 template <unsigned int N>
00101 CachingVertex<N> SequentialVertexFitter<N>::vertex(
00102     const std::vector<RefCountedVertexTrack> & tracks,
00103     const reco::BeamSpot & spot ) const
00104 {
00105   VertexState state(spot);
00106   return fit(tracks, state, true );
00107 }
00108 
00109 template <unsigned int N>
00110 CachingVertex<N> 
00111 SequentialVertexFitter<N>::vertex(const std::vector<RefCountedVertexTrack> & tracks) const
00112 {
00113   // Initial vertex state, with a very small weight matrix
00114   GlobalPoint linP = tracks[0]->linearizedTrack()->linearizationPoint();
00115   AlgebraicSymMatrix we(3,1);
00116   GlobalError error(we*10000);
00117   VertexState state(linP, error);
00118   return fit(tracks, state, false);
00119 }
00120 
00121 
00122 // Fit vertex out of a set of RecTracks. 
00123 // Uses the specified linearization point.
00124 //
00125 template <unsigned int N>
00126 CachingVertex<N>  
00127 SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack> & tracks, 
00128                                const GlobalPoint& linPoint) const
00129 { 
00130   // Initial vertex state, with a very large error matrix
00131   AlgebraicSymMatrix we(3,1);
00132   GlobalError error(we*10000);
00133   VertexState state(linPoint, error);
00134   std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
00135   return fit(vtContainer, state, false);
00136 }
00137 
00138 
00143 template <unsigned int N>
00144 CachingVertex<N> 
00145 SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack> & tracks,
00146                                const BeamSpot& beamSpot) const
00147 {
00148   VertexState beamSpotState(beamSpot);
00149   std::vector<RefCountedVertexTrack> vtContainer;
00150 
00151   if (tracks.size() > 1) {
00152     // Linearization Point search if there are more than 1 track
00153     GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
00154     if (!insideTrackerBounds(linP)) linP = GlobalPoint(0,0,0);
00155     AlgebraicSymMatrix we(3,1);
00156     GlobalError error(we*10000);
00157     VertexState lpState(linP, error);
00158     vtContainer = linearizeTracks(tracks, lpState);
00159   } else {
00160     // otherwise take the beamspot position.
00161     vtContainer = linearizeTracks(tracks, beamSpotState);
00162   }
00163 
00164   return fit(vtContainer, beamSpotState, true);
00165 }
00166 
00167 
00168 // Fit vertex out of a set of RecTracks. 
00169 // Uses the position as both the linearization point AND as prior
00170 // estimate of the vertex position. The error is used for the 
00171 // weight of the prior estimate.
00172 //
00173 template <unsigned int N>
00174 CachingVertex<N> SequentialVertexFitter<N>::vertex(
00175   const std::vector<reco::TransientTrack> & tracks, 
00176   const GlobalPoint& priorPos,
00177   const GlobalError& priorError) const
00178 { 
00179   VertexState state(priorPos, priorError);
00180   std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
00181   return fit(vtContainer, state, true);
00182 }
00183 
00184 // Fit vertex out of a set of VertexTracks
00185 // Uses the position and error for the prior estimate of the vertex.
00186 // This position is not used to relinearize the tracks.
00187 //
00188 template <unsigned int N>
00189 CachingVertex<N> SequentialVertexFitter<N>::vertex(
00190   const std::vector<RefCountedVertexTrack> & tracks, 
00191   const GlobalPoint& priorPos,
00192   const GlobalError& priorError) const
00193 {
00194   VertexState state(priorPos, priorError);
00195   return fit(tracks, state, true);
00196 }
00197 
00198 
00199 // Construct a container of VertexTrack from a set of RecTracks.
00200 //
00201 template <unsigned int N>
00202 vector<typename SequentialVertexFitter<N>::RefCountedVertexTrack> 
00203 SequentialVertexFitter<N>::linearizeTracks(
00204   const std::vector<reco::TransientTrack> & tracks, 
00205   const VertexState state) const
00206 {
00207   GlobalPoint linP = state.position();
00208   std::vector<RefCountedVertexTrack> finalTracks;
00209   finalTracks.reserve(tracks.size());
00210   for(vector<reco::TransientTrack>::const_iterator i = tracks.begin(); 
00211       i != tracks.end(); i++) {
00212     RefCountedLinearizedTrackState lTrData 
00213       = theLTrackFactory->linearizedTrackState(linP, *i);
00214     RefCountedVertexTrack vTrData = theVTrackFactory.vertexTrack(lTrData,state);
00215     finalTracks.push_back(vTrData);
00216   }
00217   return finalTracks;
00218 }
00219 
00220 
00221 // Construct new a container of VertexTrack with a new linearization point
00222 // and vertex state, from an existing set of VertexTrack, from which only the 
00223 // recTracks will be used.
00224 //
00225 template <unsigned int N>
00226 vector<typename SequentialVertexFitter<N>::RefCountedVertexTrack> 
00227 SequentialVertexFitter<N>::reLinearizeTracks(
00228   const std::vector<RefCountedVertexTrack> & tracks, 
00229   const VertexState state) const
00230 {
00231 
00232   GlobalPoint linP = state.position();
00233   std::vector<RefCountedVertexTrack> finalTracks;
00234   finalTracks.reserve(tracks.size());
00235   for(typename std::vector<RefCountedVertexTrack>::const_iterator i = tracks.begin(); 
00236       i != tracks.end(); i++) {
00237     RefCountedLinearizedTrackState lTrData = 
00238           (**i).linearizedTrack()->stateWithNewLinearizationPoint(linP);
00239     //    RefCountedLinearizedTrackState lTrData = 
00240     //      theLTrackFactory->linearizedTrackState(linP, 
00241     //                              (**i).linearizedTrack()->track());
00242     RefCountedVertexTrack vTrData =
00243       theVTrackFactory.vertexTrack(lTrData,state, (**i).weight() );
00244     finalTracks.push_back(vTrData);
00245   }
00246   return finalTracks;
00247 }
00248 
00249 
00250 // The method where the vertex fit is actually done!
00251 //
00252 template <unsigned int N>
00253 CachingVertex<N> 
00254 SequentialVertexFitter<N>::fit(const std::vector<RefCountedVertexTrack> & tracks,
00255                             const VertexState priorVertex, bool withPrior ) const
00256 {
00257   std::vector<RefCountedVertexTrack> initialTracks;
00258   GlobalPoint priorVertexPosition = priorVertex.position();
00259   GlobalError priorVertexError = priorVertex.error();
00260   
00261   CachingVertex<N> returnVertex(priorVertexPosition,priorVertexError,initialTracks,0);
00262   if (withPrior) {
00263     returnVertex = CachingVertex<N>(priorVertexPosition,priorVertexError,
00264                 priorVertexPosition,priorVertexError,initialTracks,0);
00265   }
00266   CachingVertex<N> initialVertex = returnVertex;
00267   std::vector<RefCountedVertexTrack> globalVTracks = tracks;
00268 
00269   // main loop through all the VTracks
00270   bool validVertex = true;
00271   int step = 0;
00272   GlobalPoint newPosition = priorVertexPosition;
00273   GlobalPoint previousPosition;
00274   do {
00275     CachingVertex<N> fVertex = initialVertex;
00276     // make new linearized and vertex tracks for the next iteration
00277     if(step != 0) globalVTracks = reLinearizeTracks(tracks, 
00278                                         returnVertex.vertexState());
00279 
00280     // update sequentially the vertex estimate
00281     for (typename std::vector<RefCountedVertexTrack>::const_iterator i 
00282            = globalVTracks.begin(); i != globalVTracks.end(); i++) {
00283       fVertex = theUpdator->add(fVertex,*i);
00284       if (!fVertex.isValid()) break;
00285     }
00286 
00287     validVertex = fVertex.isValid();
00288     // check tracker bounds and NaN in position
00289     if (validVertex && hasNan(fVertex.position())) {
00290       LogDebug("RecoVertex/SequentialVertexFitter") 
00291          << "Fitted position is NaN.\n";
00292       validVertex = false;
00293     }
00294 
00295     if (validVertex && !insideTrackerBounds(fVertex.position())) {
00296       LogDebug("RecoVertex/SequentialVertexFitter") 
00297          << "Fitted position is out of tracker bounds.\n";
00298       validVertex = false;
00299     }
00300 
00301     if (!validVertex) {
00302       // reset initial vertex position to (0,0,0) and force new iteration 
00303       // if number of steps not exceeded
00304       GlobalError error(AlgebraicSymMatrix(3,1)*10000);
00305       fVertex = CachingVertex<N>(GlobalPoint(0,0,0), error,
00306                               initialTracks, 0);
00307     }
00308 
00309     previousPosition = newPosition;
00310     newPosition = fVertex.position();
00311 
00312     returnVertex = fVertex;
00313     globalVTracks.clear();
00314     step++;
00315   } while ( (step != theMaxStep) &&
00316             (((previousPosition - newPosition).transverse() > theMaxShift) ||
00317                 (!validVertex) ) );
00318 
00319   if (!validVertex) {
00320     LogDebug("RecoVertex/SequentialVertexFitter") 
00321        << "Fitted position is invalid (out of tracker bounds or has NaN). Returned vertex is invalid\n";
00322     return CachingVertex<N>(); // return invalid vertex
00323   }
00324 
00325   if (step >= theMaxStep) {
00326     LogDebug("RecoVertex/SequentialVertexFitter") 
00327        << "The maximum number of steps has been exceeded. Returned vertex is invalid\n";
00328     return CachingVertex<N>(); // return invalid vertex
00329   }
00330 
00331   // smoothing
00332   returnVertex = theSmoother->smooth(returnVertex);
00333 
00334   return returnVertex;
00335 }
00336 
00337 template class SequentialVertexFitter<5>;
00338 template class SequentialVertexFitter<6>;