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