Go to the documentation of this file.00001 #include <iostream>
00002 #include <vector>
00003 #include <memory>
00004 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionHitChecker.h"
00005
00006
00007 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00008
00009
00010
00011
00012 std::pair<uint8_t,Measurement1DFloat> ConversionHitChecker::nHitsBeforeVtx(const Trajectory &traj, const reco::Vertex &vtx, double sigmaTolerance) const {
00013
00014 const std::vector<TrajectoryMeasurement> &measurements = traj.measurements();
00015
00016
00017 if (measurements.size()==0) {
00018 return std::pair<unsigned int,Measurement1DFloat>(0,Measurement1DFloat());
00019 }
00020
00021
00022 bool inout = (traj.direction() == alongMomentum);
00023
00024 GlobalPoint vtxPos(vtx.x(),vtx.y(),vtx.z());
00025
00026 uint8_t nhits = 0;
00027
00028 std::vector<TrajectoryMeasurement>::const_iterator closest = measurements.begin();
00029 double distance = 1e6;
00030
00031 for (std::vector<TrajectoryMeasurement>::const_iterator fwdit = measurements.begin(), revit = measurements.end()-1;
00032 fwdit != measurements.end(); ++fwdit,--revit) {
00033
00034 std::vector<TrajectoryMeasurement>::const_iterator it;
00035 if (inout) {
00036 it = fwdit;
00037 }
00038 else {
00039 it = revit;
00040 }
00041
00042
00043 if (it->recHit() && it->recHit()->isValid()) {
00044 ++nhits;
00045 distance = (vtxPos - it->updatedState().globalPosition()).mag();
00046 closest = it;
00047 }
00048
00049 if ( (measurements.end()-fwdit)==1) {
00050 break;
00051 }
00052
00053
00054 std::vector<TrajectoryMeasurement>::const_iterator nextit;
00055 if (inout) {
00056 nextit = it+1;
00057 }
00058 else {
00059 nextit = it-1;
00060 }
00061
00062
00063 if ( nextit->recHit() && nextit->recHit()->isValid() ) {
00064 double nextDistance = (vtxPos - nextit->updatedState().globalPosition()).mag();
00065 if (nextDistance > distance) {
00066 break;
00067 }
00068 }
00069
00070 }
00071
00072
00073
00074
00075 GlobalVector momDir = closest->updatedState().globalMomentum().unit();
00076 double decayLengthHitToVtx = (vtxPos - closest->updatedState().globalPosition()).dot(momDir);
00077
00078 AlgebraicVector3 j;
00079 j[0] = momDir.x();
00080 j[1] = momDir.y();
00081 j[2] = momDir.z();
00082 AlgebraicVector6 jj;
00083 jj[0] = momDir.x();
00084 jj[1] = momDir.y();
00085 jj[2] = momDir.z();
00086 jj[3] =0.;
00087 jj[4] =0.;
00088 jj[5] =0.;
00089
00090
00091
00092 double trackError2 = ROOT::Math::Similarity(jj,closest->updatedState().cartesianError().matrix());
00093 double vertexError2 = ROOT::Math::Similarity(j,vtx.covariance());
00094 double decayLenError = sqrt(trackError2+vertexError2);
00095
00096 Measurement1DFloat decayLength(decayLengthHitToVtx,decayLenError);
00097
00098 if (decayLength.significance() < sigmaTolerance) {
00099
00100 --nhits;
00101 }
00102
00103 return std::pair<unsigned int,Measurement1DFloat>(nhits,decayLength);
00104
00105 }
00106
00107 uint8_t ConversionHitChecker::nSharedHits(const reco::Track &trk1, const reco::Track &trk2) const {
00108
00109 uint8_t nShared = 0;
00110
00111 for (trackingRecHit_iterator iHit1 = trk1.recHitsBegin(); iHit1 != trk1.recHitsEnd(); ++iHit1) {
00112 const TrackingRecHit *hit1 = iHit1->get();
00113 if (hit1->isValid()) {
00114 for (trackingRecHit_iterator iHit2 = trk2.recHitsBegin(); iHit2 != trk2.recHitsEnd(); ++iHit2) {
00115 const TrackingRecHit *hit2 = iHit2->get();
00116 if (hit2->isValid() && hit1->sharesInput(hit2,TrackingRecHit::some)) {
00117 ++nShared;
00118 }
00119 }
00120 }
00121 }
00122
00123 return nShared;
00124
00125 }
00126
00127