CMS 3D CMS Logo

VertexFitterResult Class Reference

Very basic class containing only the positions of the simulated and reconstructed vertices, total chi**2, chi**2 probability and number of degrees of freedom. More...

#include <RecoVertex/KalmanVertexFit/interface/VertexFitterResult.h>

List of all members.

Public Types

typedef std::vector
< reco::TransientTrack
TTrackCont

Public Member Functions

const float * chi2Information () const
void fill (const TransientVertex &recVertex, const TTrackCont &recTrackV, const TrackingVertex *simv=0, reco::RecoToSimCollection *recSimColl=0, const float &time=0)
void fill (const TransientVertex &recv, const TrackingVertex *simv=0, reco::RecoToSimCollection *recSimColl=0, const float &time=0)
const intnumberRecTracks ()
const intnumberSimTracks ()
const float * recErrors (const int i) const
const float * recParameters (const int i) const
const intrecTrack_simIndex ()
const float * recTrackWeight ()
const float * recVertexErr () const
const float * recVertexPos () const
const float * refErrors (const int i) const
const float * refParameters (const int i) const
void reset ()
const float * simParameters (const int i) const
const intsimTrack_recIndex ()
const float * simVertexPos () const
const float * time () const
const inttrackInformation () const
 VertexFitterResult (const int maxTracks=100, TrackAssociatorByChi2 *associator=0)
const intvertexPresent () const
 ~VertexFitterResult ()

Private Member Functions

void fillErrors (const PerigeeTrajectoryError &pte, float *errors[5], int trackNumber)
void fillParameters (const PerigeeTrajectoryParameters &ptp, float *params[5], int trackNumber)
void fillParameters (const reco::TrackBase::ParameterVector &params, float *params[5], int trackNumber)

Private Attributes

TrackAssociatorByChi2associatorForParamAtPca
float chi [3]
float fitTime
int numberOfRecTracks
int numberOfSimTracks
float recErr [3]
float * recErrs [5]
intrecIndex
float * recPars [5]
float recPos [3]
float * refErrs [5]
float * refPars [5]
intsimIndex
float * simPars [5]
float simPos [3]
int theMaxTracks
int tracks [3]
float * trackWeight
int vertex


Detailed Description

Very basic class containing only the positions of the simulated and reconstructed vertices, total chi**2, chi**2 probability and number of degrees of freedom.

The only thing to be done is to call the method fill for each vertex.

Definition at line 24 of file VertexFitterResult.h.


Member Typedef Documentation

typedef std::vector<reco::TransientTrack> VertexFitterResult::TTrackCont

Definition at line 28 of file VertexFitterResult.h.


Constructor & Destructor Documentation

VertexFitterResult::VertexFitterResult ( const int  maxTracks = 100,
TrackAssociatorByChi2 associator = 0 
)

Definition at line 9 of file VertexFitterResult.cc.

References associatorForParamAtPca, i, numberOfRecTracks, numberOfSimTracks, recErrs, recIndex, recPars, refErrs, refPars, reset(), simIndex, simPars, theMaxTracks, and trackWeight.

00010         : associatorForParamAtPca(associator)
00011 {
00012   theMaxTracks = maxTracks;
00013   if (associatorForParamAtPca==0) theMaxTracks=0;
00014   for ( int i=0; i<5; i++ ) {
00015     if ( maxTracks>0 ) {
00016       simPars[i] = new float[maxTracks];
00017       recPars[i] = new float[maxTracks];
00018       refPars[i] = new float[maxTracks];
00019       recErrs[i] = new float[maxTracks];
00020       refErrs[i] = new float[maxTracks];
00021     }
00022     else {
00023       simPars[i] = 0;
00024       recPars[i] = 0;
00025       refPars[i] = 0;
00026       recErrs[i] = 0;
00027       refErrs[i] = 0;
00028     }
00029   }
00030   trackWeight = new float[maxTracks];
00031   simIndex = new int[maxTracks];
00032   recIndex = new int[maxTracks];
00033   numberOfRecTracks=theMaxTracks;
00034   numberOfSimTracks=theMaxTracks;
00035   reset();
00036 }

VertexFitterResult::~VertexFitterResult (  ) 

Definition at line 38 of file VertexFitterResult.cc.

References i, recErrs, recIndex, recPars, refErrs, refPars, simIndex, simPars, and trackWeight.

00039 {
00040     //
00041     // delete arrays
00042     //
00043     for ( int i=0; i<5; i++ ) {
00044       delete [] simPars[i];
00045       delete [] recPars[i];
00046       delete [] refPars[i];
00047       delete [] recErrs[i];
00048       delete [] refErrs[i];
00049     }
00050     delete trackWeight;
00051     delete simIndex;
00052     delete recIndex;
00053 }


Member Function Documentation

const float* VertexFitterResult::chi2Information (  )  const [inline]

Definition at line 46 of file VertexFitterResult.h.

References chi.

Referenced by SimpleVertexTree::SimpleVertexTree().

00046 {return chi;}

void VertexFitterResult::fill ( const TransientVertex recVertex,
const TTrackCont recTrackV,
const TrackingVertex simv = 0,
reco::RecoToSimCollection recSimColl = 0,
const float &  time = 0 
)

Definition at line 64 of file VertexFitterResult.cc.

References associatorForParamAtPca, edm::RefVector< C, T, F >::begin(), chi, ChiSquaredProbability(), GlobalErrorBase< T, ErrorWeightType >::cxx(), GlobalErrorBase< T, ErrorWeightType >::cyy(), GlobalErrorBase< T, ErrorWeightType >::czz(), TrackingVertex::daughterTracks(), TrackingVertex::daughterTracks_begin(), TrackingVertex::daughterTracks_end(), TransientVertex::degreesOfFreedom(), e, edm::RefVector< C, T, F >::end(), fillErrors(), fillParameters(), find(), fitTime, TransientVertex::isValid(), numberOfRecTracks, numberOfSimTracks, TransientVertex::originalTracks(), TrackAssociatorByChi2::parametersAtClosestApproach(), TrajectoryStateClosestToPoint::perigeeError(), TrackingVertex::position(), TransientVertex::position(), TransientVertex::positionError(), recErr, recErrs, recIndex, recPars, recPos, simIndex, simPars, simPos, funct::sqrt(), theMaxTracks, TransientVertex::totalChiSquared(), reco::TrackTransientTrack::trackBaseRef(), tracks, vertex, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

00067 {
00068   TrackingParticleRefVector simTrackV;
00069 
00070   Basic3DVector<double> vert;
00071   if (recVertex.isValid()) {
00072     recPos[0] = recVertex.position().x();
00073     recPos[1] = recVertex.position().y();
00074     recPos[2] = recVertex.position().z();
00075 
00076     recErr[0] = sqrt(recVertex.positionError().cxx());
00077     recErr[1] = sqrt(recVertex.positionError().cyy());
00078     recErr[2] = sqrt(recVertex.positionError().czz());
00079     vert = (Basic3DVector<double>) recVertex.position();
00080 
00081     chi[0] = recVertex.totalChiSquared();
00082     chi[1] = recVertex.degreesOfFreedom();
00083     chi[2] = ChiSquaredProbability(recVertex.totalChiSquared(), 
00084                                            recVertex.degreesOfFreedom());
00085     vertex = 2;
00086     fitTime = time;
00087     tracks[1] = recVertex.originalTracks().size();
00088   }
00089 
00090   if (simv!=0) {
00091     simPos[0] = simv->position().x();
00092     simPos[1] = simv->position().y();
00093     simPos[2] = simv->position().z();
00094 
00095     simTrackV = simv->daughterTracks();
00096     vertex += 1;
00097     for(TrackingVertex::tp_iterator simTrack = simv->daughterTracks_begin();
00098                  (simTrack != simv->daughterTracks_end() && (numberOfSimTracks<theMaxTracks));
00099                  simTrack++) {
00100       
00101       Basic3DVector<double> momAtVtx((**simTrack).momentum());
00102 
00103       pair<bool, reco::TrackBase::ParameterVector> paramPair =
00104         associatorForParamAtPca->parametersAtClosestApproach(vert, momAtVtx, 
00105           (float) (**simTrack).charge(), recTrackV.front().stateAtBeamLine().beamSpot());
00106         if (paramPair.first) {
00107           fillParameters(paramPair.second, simPars, numberOfSimTracks);
00108           simIndex[numberOfSimTracks] = -1;
00109           ++numberOfSimTracks;
00110         }
00111     }
00112     tracks[0] = numberOfSimTracks;
00113   }
00114 
00115 
00116   // now store all the recTrack...  
00117 
00118   for(TTrackCont::const_iterator recTrack =recTrackV.begin();
00119                (recTrack != recTrackV.end() 
00120                 && (numberOfRecTracks<theMaxTracks));
00121                recTrack++) {
00122     //    cout << "Input; 1/Pt " << 1./(*recTrack).momentumAtVertex().transverse() << endl;
00123 
00124     //looking for sim tracks corresponding to our reconstructed tracks:  
00125     simIndex[numberOfRecTracks] = -1;
00126 
00127     std::vector<std::pair<TrackingParticleRef, double> > simFound;
00128     try {
00129       const TrackTransientTrack* ttt = dynamic_cast<const TrackTransientTrack*>(recTrack->basicTransientTrack());
00130       if ((ttt!=0) && (recSimColl!=0)) simFound = (*recSimColl)[ttt->trackBaseRef()];
00131 //       if (recSimColl!=0) simFound = (*recSimColl)[recTrack->persistentTrackRef()];
00132 //      if (recSimColl!=0) simFound = (*recSimColl)[recTrack];
00133 
00134     } catch (cms::Exception e) {
00135 //       LogDebug("TrackValidator") << "reco::Track #" << rT << " with pt=" << track->pt() 
00136 //                               << " NOT associated to any TrackingParticle" << "\n";
00137 //       edm::LogError("TrackValidator") << e.what() << "\n";
00138     }
00139 
00140     if(simFound.size() != 0) {
00141       //OK, it was associated, so get the state on the same surface as the 'SimState'
00142       TrackingParticleRefVector::const_iterator simTrackI = 
00143         find(simTrackV.begin(), simTrackV.end(), simFound[0].first);
00144       if (simTrackI!=simTrackV.end()) ++tracks[2];
00145       int simTrackIndex = simTrackI-simTrackV.begin();
00146       if (simTrackIndex<numberOfSimTracks) {
00147         simIndex[numberOfRecTracks] = simTrackIndex;
00148         recIndex[simTrackIndex] = numberOfRecTracks;
00149         //      cout << "Assoc; 1/Pt " << 1./(*recTrack).momentumAtVertex().transverse() << endl;
00150       }
00151     }
00152 
00153     TrajectoryStateClosestToPoint tscp = recTrack->trajectoryStateClosestToPoint(recVertex.position());
00154     fillParameters(recTrack->track().parameters(), recPars, numberOfRecTracks);
00155     fillErrors(tscp.perigeeError(), recErrs, numberOfRecTracks);
00156 //     trackWeight[numberOfRecTracks] = recVertex.trackWeight(*recTrack);
00157 // 
00158 //     if ((recVertex.isValid())&&(recVertex.hasRefittedTracks())) {
00159 //       //looking for corresponding refitted tracks:
00160 //       TrajectoryStateOnSurface refip;
00161 //       RecTrackCont::iterator refTrackI = 
00162 //                      find_if(refTracks.begin(), refTracks.end(), RecTrackMatch(*recTrack));
00163 //       if (refTrackI!=refTracks.end()) {
00164 //         // If it was not found, it would mean that it was not used in the fit,
00165 //      // or with a low weight such that the track was then discarded.
00166 //      if(simFound.size() != 0) {
00167 //        refip = refTrackI->stateOnSurface(simFound[0]->impactPointStateOnSurface().surface());
00168 //      } else {
00169 //           refip = refTrackI->innermostState();
00170 //      }
00171 // 
00172 //      fillParameters(refip, refPars, numberOfRecTracks);
00173 //      fillErrors(refip, refErrs, numberOfRecTracks);
00174 //       }
00175 //     }
00176 // 
00177     ++numberOfRecTracks;
00178   }
00179   
00180 }

void VertexFitterResult::fill ( const TransientVertex recv,
const TrackingVertex simv = 0,
reco::RecoToSimCollection recSimColl = 0,
const float &  time = 0 
)

Definition at line 55 of file VertexFitterResult.cc.

References TransientVertex::isValid(), and TransientVertex::originalTracks().

Referenced by SimpleVertexTree::fill().

00058 {
00059   TTrackCont recTrackV;
00060   if (recVertex.isValid()) recTrackV = recVertex.originalTracks();
00061   fill(recVertex, recTrackV, simv, recSimColl, time);
00062 }

void VertexFitterResult::fillErrors ( const PerigeeTrajectoryError pte,
float *  errors[5],
int  trackNumber 
) [private]

Definition at line 203 of file VertexFitterResult.cc.

References PerigeeTrajectoryError::longitudinalImpactParameterError(), PerigeeTrajectoryError::phiError(), PerigeeTrajectoryError::thetaError(), PerigeeTrajectoryError::transverseCurvatureError(), and PerigeeTrajectoryError::transverseImpactParameterError().

Referenced by fill().

00205 {
00206   errors[0][trackNumber] = pte.transverseCurvatureError(); 
00207   errors[1][trackNumber] = pte.thetaError(); 
00208   errors[2][trackNumber] = pte.phiError(); 
00209   errors[3][trackNumber] = pte.transverseImpactParameterError(); 
00210   errors[4][trackNumber] = pte.longitudinalImpactParameterError(); 
00211 }

void VertexFitterResult::fillParameters ( const PerigeeTrajectoryParameters ptp,
float *  params[5],
int  trackNumber 
) [private]

Definition at line 192 of file VertexFitterResult.cc.

References PerigeeTrajectoryParameters::vector_old().

00194 {
00195   const AlgebraicVector & perigee = ptp.vector_old();
00196   params[0][trackNumber] = perigee[0];
00197   params[1][trackNumber] = perigee[1];
00198   params[2][trackNumber] = perigee[2];
00199   params[3][trackNumber] = perigee[3];
00200   params[4][trackNumber] = perigee[4];
00201 }

void VertexFitterResult::fillParameters ( const reco::TrackBase::ParameterVector params,
float *  params[5],
int  trackNumber 
) [private]

Definition at line 182 of file VertexFitterResult.cc.

Referenced by fill().

00184 {
00185   params[0][trackNumber] = perigee[0];
00186   params[1][trackNumber] = perigee[1];
00187   params[2][trackNumber] = perigee[2];
00188   params[3][trackNumber] = perigee[3];
00189   params[4][trackNumber] = perigee[4];
00190 }

const int* VertexFitterResult::numberRecTracks (  )  [inline]

Definition at line 51 of file VertexFitterResult.h.

References numberOfRecTracks.

Referenced by SimpleVertexTree::SimpleVertexTree().

00051 {return &numberOfRecTracks;}

const int* VertexFitterResult::numberSimTracks (  )  [inline]

Definition at line 50 of file VertexFitterResult.h.

References numberOfSimTracks.

Referenced by SimpleVertexTree::SimpleVertexTree().

00050 {return &numberOfSimTracks;}

const float* VertexFitterResult::recErrors ( const int  i  )  const [inline]

Definition at line 72 of file VertexFitterResult.h.

References recErrs.

Referenced by SimpleVertexTree::SimpleVertexTree().

00073   {
00074     if ( i<0 || i>=5 )  return 0;
00075     return recErrs[i];
00076   }

const float* VertexFitterResult::recParameters ( const int  i  )  const [inline]

Definition at line 57 of file VertexFitterResult.h.

References recPars.

Referenced by SimpleVertexTree::SimpleVertexTree().

00058   {
00059     if ( i<0 || i>=5 )  return 0;
00060     return recPars[i];
00061   }

const int* VertexFitterResult::recTrack_simIndex (  )  [inline]

Definition at line 55 of file VertexFitterResult.h.

References simIndex.

Referenced by SimpleVertexTree::SimpleVertexTree().

00055 {return simIndex;}

const float* VertexFitterResult::recTrackWeight (  )  [inline]

Definition at line 56 of file VertexFitterResult.h.

References trackWeight.

Referenced by SimpleVertexTree::SimpleVertexTree().

00056 {return trackWeight;}

const float* VertexFitterResult::recVertexErr (  )  const [inline]

Definition at line 44 of file VertexFitterResult.h.

References recErr.

Referenced by SimpleVertexTree::SimpleVertexTree().

00044 {return recErr;}

const float* VertexFitterResult::recVertexPos (  )  const [inline]

Definition at line 43 of file VertexFitterResult.h.

References recPos.

Referenced by SimpleVertexTree::SimpleVertexTree().

00043 {return recPos;}

const float* VertexFitterResult::refErrors ( const int  i  )  const [inline]

Definition at line 77 of file VertexFitterResult.h.

References refErrs.

Referenced by SimpleVertexTree::SimpleVertexTree().

00078   {
00079     if ( i<0 || i>=5 )  return 0;
00080     return refErrs[i];
00081   }

const float* VertexFitterResult::refParameters ( const int  i  )  const [inline]

Definition at line 62 of file VertexFitterResult.h.

References refPars.

Referenced by SimpleVertexTree::SimpleVertexTree().

00063   {
00064     if ( i<0 || i>=5 )  return 0;
00065     return refPars[i];
00066   }

void VertexFitterResult::reset ( void   ) 

Definition at line 213 of file VertexFitterResult.cc.

References chi, fitTime, i, j, numberOfRecTracks, numberOfSimTracks, recErr, recErrs, recIndex, recPars, recPos, refErrs, refPars, simIndex, simPars, simPos, tracks, trackWeight, and vertex.

Referenced by SimpleVertexTree::fill(), and VertexFitterResult().

00214 {
00215   for ( int i=0; i<3; ++i ) {
00216     simPos[i] = 0.;
00217     recPos[i] = 0.;
00218     recErr[i] = 0.;
00219     chi[i] = 0.;
00220     tracks[i] = 0;
00221   }
00222   vertex =0;
00223   fitTime = 0;
00224 
00225   for ( int j=0; j<numberOfRecTracks; ++j ) {
00226     for ( int i=0; i<5; ++i ) {
00227        recPars[i][j] = 0;
00228        refPars[i][j] = 0;
00229        recErrs[i][j] = 0;
00230        refErrs[i][j] = 0;
00231     }
00232     trackWeight[j] = 0;
00233     simIndex[j] = -1;
00234   }
00235   for ( int j=0; j<numberOfSimTracks; ++j ) {
00236     for ( int i=0; i<5; ++i ) {
00237        simPars[i][j] = 0;
00238     }
00239     recIndex[j] = -1;
00240   }
00241 
00242   numberOfRecTracks=0;
00243   numberOfSimTracks=0;
00244 }

const float* VertexFitterResult::simParameters ( const int  i  )  const [inline]

Definition at line 67 of file VertexFitterResult.h.

References simPars.

Referenced by SimpleVertexTree::SimpleVertexTree().

00068   {
00069     if ( i<0 || i>=5 )  return 0;
00070     return simPars[i];
00071   }

const int* VertexFitterResult::simTrack_recIndex (  )  [inline]

Definition at line 53 of file VertexFitterResult.h.

References recIndex.

Referenced by SimpleVertexTree::SimpleVertexTree().

00053 {return recIndex;}

const float* VertexFitterResult::simVertexPos (  )  const [inline]

Definition at line 42 of file VertexFitterResult.h.

References simPos.

Referenced by SimpleVertexTree::SimpleVertexTree().

00042 {return simPos;}

const float* VertexFitterResult::time ( void   )  const [inline]

Definition at line 48 of file VertexFitterResult.h.

References fitTime.

Referenced by SimpleVertexTree::SimpleVertexTree().

00048 {return &fitTime;}

const int* VertexFitterResult::trackInformation (  )  const [inline]

Definition at line 45 of file VertexFitterResult.h.

References tracks.

Referenced by SimpleVertexTree::SimpleVertexTree().

00045 {return tracks;}

const int* VertexFitterResult::vertexPresent (  )  const [inline]

Definition at line 47 of file VertexFitterResult.h.

References vertex.

Referenced by SimpleVertexTree::SimpleVertexTree().

00047 {return &vertex;}


Member Data Documentation

TrackAssociatorByChi2* VertexFitterResult::associatorForParamAtPca [private]

Definition at line 114 of file VertexFitterResult.h.

Referenced by fill(), and VertexFitterResult().

float VertexFitterResult::chi[3] [private]

Definition at line 119 of file VertexFitterResult.h.

Referenced by chi2Information(), fill(), and reset().

float VertexFitterResult::fitTime [private]

Definition at line 122 of file VertexFitterResult.h.

Referenced by fill(), reset(), and time().

int VertexFitterResult::numberOfRecTracks [private]

Definition at line 130 of file VertexFitterResult.h.

Referenced by fill(), numberRecTracks(), reset(), and VertexFitterResult().

int VertexFitterResult::numberOfSimTracks [private]

Definition at line 130 of file VertexFitterResult.h.

Referenced by fill(), numberSimTracks(), reset(), and VertexFitterResult().

float VertexFitterResult::recErr[3] [private]

Definition at line 118 of file VertexFitterResult.h.

Referenced by fill(), recVertexErr(), and reset().

float* VertexFitterResult::recErrs[5] [private]

Definition at line 128 of file VertexFitterResult.h.

Referenced by fill(), recErrors(), reset(), VertexFitterResult(), and ~VertexFitterResult().

int * VertexFitterResult::recIndex [private]

Definition at line 132 of file VertexFitterResult.h.

Referenced by fill(), reset(), simTrack_recIndex(), VertexFitterResult(), and ~VertexFitterResult().

float* VertexFitterResult::recPars[5] [private]

Definition at line 126 of file VertexFitterResult.h.

Referenced by fill(), recParameters(), reset(), VertexFitterResult(), and ~VertexFitterResult().

float VertexFitterResult::recPos[3] [private]

Definition at line 117 of file VertexFitterResult.h.

Referenced by fill(), recVertexPos(), and reset().

float* VertexFitterResult::refErrs[5] [private]

Definition at line 129 of file VertexFitterResult.h.

Referenced by refErrors(), reset(), VertexFitterResult(), and ~VertexFitterResult().

float* VertexFitterResult::refPars[5] [private]

Definition at line 127 of file VertexFitterResult.h.

Referenced by refParameters(), reset(), VertexFitterResult(), and ~VertexFitterResult().

int* VertexFitterResult::simIndex [private]

Definition at line 132 of file VertexFitterResult.h.

Referenced by fill(), recTrack_simIndex(), reset(), VertexFitterResult(), and ~VertexFitterResult().

float* VertexFitterResult::simPars[5] [private]

Definition at line 125 of file VertexFitterResult.h.

Referenced by fill(), reset(), simParameters(), VertexFitterResult(), and ~VertexFitterResult().

float VertexFitterResult::simPos[3] [private]

Definition at line 116 of file VertexFitterResult.h.

Referenced by fill(), reset(), and simVertexPos().

int VertexFitterResult::theMaxTracks [private]

Definition at line 124 of file VertexFitterResult.h.

Referenced by fill(), and VertexFitterResult().

int VertexFitterResult::tracks[3] [private]

Definition at line 120 of file VertexFitterResult.h.

Referenced by fill(), reset(), and trackInformation().

float* VertexFitterResult::trackWeight [private]

Definition at line 131 of file VertexFitterResult.h.

Referenced by recTrackWeight(), reset(), VertexFitterResult(), and ~VertexFitterResult().

int VertexFitterResult::vertex [private]

Definition at line 121 of file VertexFitterResult.h.

Referenced by fill(), reset(), and vertexPresent().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:35:00 2009 for CMSSW by  doxygen 1.5.4