CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoEgamma/EgammaTools/interface/ConversionFinder.h

Go to the documentation of this file.
00001 #ifndef EgammaTools_ConversionFinder_h
00002 #define EgammaTools_ConversionFinder_h
00003 
00014 #include "DataFormats/TrackReco/interface/Track.h"
00015 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00016 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00017 #include "DataFormats/Math/interface/Point3D.h"
00018 #include "MagneticField/Engine/interface/MagneticField.h"
00019 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00020 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00021 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00022 #include "DataFormats/Math/interface/LorentzVector.h"
00023 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00024 #include "Math/VectorUtil.h"
00025 #include "ConversionInfo.h"
00026 
00027 
00028 
00029 /*
00030    Class Looks for oppositely charged track in the
00031    track collection with the minimum delta cot(theta) between the track
00032    and the electron's CTF track (if it doesn't exist, we use the electron's
00033    GSF track). Calculate the dist, dcot, point of conversion and the
00034    radius of conversion for this pair and fill the ConversionInfo
00035 */
00036 
00037 class ConversionFinder {
00038  public:
00039   ConversionFinder();
00040   ~ConversionFinder();
00041   //bField has to be supplied in Tesla
00042 
00043   //returns a vector of Conversion Infos.
00044   //you have to decide which is the best
00045   std::vector<ConversionInfo> getConversionInfos(const reco::GsfElectronCore&,
00046                                                 const edm::Handle<reco::TrackCollection>& ctftracks_h,
00047                                                 const edm::Handle<reco::GsfTrackCollection>& gsftracks_h,
00048                                                 const double bFieldAtOrigin,
00049                                                 const double minFracSharedHits = 0.45);
00050 
00051   //retruns the "best" Conversion Info after calling getConversionInfos
00052   ConversionInfo getConversionInfo(const reco::GsfElectronCore&,
00053                                    const edm::Handle<reco::TrackCollection>& ctftracks_h,
00054                                    const edm::Handle<reco::GsfTrackCollection>& gsftracks_h,
00055                                    const double bFieldAtOrigin,
00056                                    const double minFracSharedHits = 0.45);
00057 
00058   //retruns the "best" Conversion Info after calling getConversionInfos
00059   ConversionInfo getConversionInfo(const reco::GsfElectron& gsfElectron,
00060                                    const edm::Handle<reco::TrackCollection>& ctftracks_h,
00061                                    const edm::Handle<reco::GsfTrackCollection>& gsftracks_h,
00062                                    const double bFieldAtOrigin,
00063                                    const double minFracSharedHits = 0.45);
00064 
00065   //used internally, so call this only if you know what you're doing.
00066   ConversionInfo getConversionInfo(const reco::Track *el_track,
00067                                    const reco::Track *candPartnerTk,
00068                                    const double bFieldAtOrigin);
00069 
00070   const reco::Track* getElectronTrack(const reco::GsfElectron &, const float minFracSharedHits = 0.45);
00071 
00072   const reco::Track* getElectronTrack(const reco::GsfElectronCore &, const float minFracSharedHits = 0.45);
00073 
00074 
00075   //takes in a vector of candidate conversion partners
00076   //and arbitrates between them returning the one with the
00077   //smallest R=sqrt(dist*dist + dcot*dcot)
00078   ConversionInfo arbitrateConversionPartnersbyR(const std::vector<ConversionInfo>& v_convCandidates);
00079 
00080 
00081   //places different cuts on dist, dcot, delmissing hits and arbitration based on R = sqrt(dist*dist + dcot*dcot)
00082   ConversionInfo findBestConversionMatch(const std::vector<ConversionInfo>& v_convCandidates);
00083 
00084   //function below is only for backwards compatibility
00085   static std::pair<double, double> getConversionInfo(math::XYZTLorentzVector trk1_p4,
00086                                                      int trk1_q, float trk1_d0,
00087                                                      math::XYZTLorentzVector trk2_p4,
00088                                                      int trk2_q, float trk2_d0,
00089                                                      float bFieldAtOrigin);
00090 
00091   //for backwards compatibility. Does not use the GSF track collection. This function will be
00092   //deprecated soon
00093   ConversionInfo getConversionInfo(const reco::GsfElectron& gsfElectron,
00094                                    const edm::Handle<reco::TrackCollection>& track_h,
00095                                    const double bFieldAtOrigin,
00096                                    const double minFracSharedHits = 0.45);
00097 
00098 
00099   // DEPRECATED
00100   bool isFromConversion( const ConversionInfo &, double maxAbsDist = 0.02, double maxAbsDcot = 0.02);
00101 
00102  private:
00103 
00104 };
00105 #endif