CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoVertex/KalmanVertexFit/src/VertexFitterResult.cc

Go to the documentation of this file.
00001 #include "RecoVertex/KalmanVertexFit/interface/VertexFitterResult.h"
00002 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00003 #include "SimTracker/TrackAssociation/interface/TrackAssociatorBase.h"
00004 #include "TrackingTools/TransientTrack/interface/TrackTransientTrack.h"
00005 
00006 using namespace reco;
00007 using namespace std;
00008 
00009 VertexFitterResult::VertexFitterResult(const int maxTracks, TrackAssociatorByChi2 *associator)
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 }
00037 
00038 VertexFitterResult::~VertexFitterResult()
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 }
00054 
00055 void VertexFitterResult::fill(const TransientVertex & recVertex,
00056         const TrackingVertex * simv, reco::RecoToSimCollection *recSimColl,
00057         const float &time) 
00058 {
00059   TTrackCont recTrackV;
00060   if (recVertex.isValid()) recTrackV = recVertex.originalTracks();
00061   fill(recVertex, recTrackV, simv, recSimColl, time);
00062 }
00063 
00064 void VertexFitterResult::fill(const TransientVertex & recVertex, 
00065         const TTrackCont & recTrackV, const TrackingVertex * simv,
00066         reco::RecoToSimCollection *recSimColl, const float &time)
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       std::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     //    std::cout << "Input; 1/Pt " << 1./(*recTrack).momentumAtVertex().transverse() << std::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() << std::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 }
00181 
00182 void VertexFitterResult::fillParameters (const reco::TrackBase::ParameterVector& perigee,
00183         float* params[5], int trackNumber)
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 }
00191 
00192 void VertexFitterResult::fillParameters (const PerigeeTrajectoryParameters & ptp,
00193         float* params[5], int trackNumber)
00194 {
00195   const AlgebraicVector5 & perigee = ptp.vector();
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 }
00202 
00203 void VertexFitterResult::fillErrors (const PerigeeTrajectoryError & pte,
00204         float* errors[5], int trackNumber)
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 }
00212 
00213 void VertexFitterResult::reset()
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 }