Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:42

0001 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0002 #include "DataFormats/Common/interface/Handle.h"
0003 #include "TrackingTools/TransientTrack/interface/GsfTransientTrack.h"
0004 #include "TrackingTools/TransientTrack/interface/TransientTrackFromFTS.h"
0005 #include "FWCore/Utilities/interface/isFinite.h"
0006 
0007 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0008 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
0009 
0010 #include "DataFormats/PatCandidates/interface/PackedCandidate.h"
0011 
0012 using namespace reco;
0013 using namespace std;
0014 using namespace edm;
0015 
0016 TransientTrack TransientTrackBuilder::build(const Track* t) const {
0017   return TransientTrack(*t, theField, theTrackingGeometry);
0018 }
0019 
0020 TransientTrack TransientTrackBuilder::build(const Track& t) const {
0021   return TransientTrack(t, theField, theTrackingGeometry);
0022 }
0023 
0024 TransientTrack TransientTrackBuilder::build(const GsfTrack* t) const {
0025   return TransientTrack(new GsfTransientTrack(*t, theField, theTrackingGeometry));
0026 }
0027 
0028 TransientTrack TransientTrackBuilder::build(const GsfTrack& t) const {
0029   return TransientTrack(new GsfTransientTrack(t, theField, theTrackingGeometry));
0030 }
0031 
0032 TransientTrack TransientTrackBuilder::build(const CandidatePtr* t) const {
0033   reco::PFCandidatePtr tryPF(*t);
0034   edm::Ptr<pat::PackedCandidate> tryPacked(*t);
0035   if (tryPF.get() != nullptr && tryPF->isTimeValid()) {
0036     return TransientTrack(*t, tryPF->time(), tryPF->timeError(), theField, theTrackingGeometry);
0037   } else if (tryPacked.get() != nullptr && tryPacked->timeError() > 0.f) {
0038     return TransientTrack(*t, (double)tryPacked->time(), (double)tryPacked->timeError(), theField, theTrackingGeometry);
0039   }
0040   return TransientTrack(*t, theField, theTrackingGeometry);
0041 }
0042 
0043 TransientTrack TransientTrackBuilder::build(const CandidatePtr& t) const { return this->build(&t); }
0044 
0045 TransientTrack TransientTrackBuilder::build(const TrackRef* t) const {
0046   return TransientTrack(*t, theField, theTrackingGeometry);
0047 }
0048 
0049 TransientTrack TransientTrackBuilder::build(const TrackRef& t) const {
0050   return TransientTrack(t, theField, theTrackingGeometry);
0051 }
0052 
0053 TransientTrack TransientTrackBuilder::build(const GsfTrackRef* t) const {
0054   return TransientTrack(new GsfTransientTrack(*t, theField, theTrackingGeometry));
0055 }
0056 
0057 TransientTrack TransientTrackBuilder::build(const GsfTrackRef& t) const {
0058   return TransientTrack(new GsfTransientTrack(t, theField, theTrackingGeometry));
0059 }
0060 
0061 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<reco::TrackCollection>& trkColl) const {
0062   vector<TransientTrack> ttVect;
0063   ttVect.reserve((*trkColl).size());
0064   for (unsigned int i = 0; i < (*trkColl).size(); i++) {
0065     ttVect.push_back(TransientTrack(TrackRef(trkColl, i), theField, theTrackingGeometry));
0066   }
0067   return ttVect;
0068 }
0069 
0070 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<reco::GsfTrackCollection>& trkColl) const {
0071   vector<TransientTrack> ttVect;
0072   ttVect.reserve((*trkColl).size());
0073   for (unsigned int i = 0; i < (*trkColl).size(); i++) {
0074     ttVect.push_back(TransientTrack(new GsfTransientTrack(GsfTrackRef(trkColl, i), theField, theTrackingGeometry)));
0075   }
0076   return ttVect;
0077 }
0078 
0079 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<edm::View<Track> >& trkColl) const {
0080   vector<TransientTrack> ttVect;
0081   ttVect.reserve((*trkColl).size());
0082   for (unsigned int i = 0; i < (*trkColl).size(); i++) {
0083     const Track* trk = &(*trkColl)[i];
0084     const GsfTrack* gsfTrack = dynamic_cast<const GsfTrack*>(trk);
0085     if (gsfTrack) {
0086       ttVect.push_back(TransientTrack(
0087           new GsfTransientTrack(RefToBase<Track>(trkColl, i).castTo<GsfTrackRef>(), theField, theTrackingGeometry)));
0088     } else {  // gsf
0089       ttVect.push_back(TransientTrack(RefToBase<Track>(trkColl, i).castTo<TrackRef>(), theField, theTrackingGeometry));
0090     }
0091   }
0092   return ttVect;
0093 }
0094 
0095 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<reco::TrackCollection>& trkColl,
0096                                                     const edm::ValueMap<float>& trackTimes,
0097                                                     const edm::ValueMap<float>& trackTimeResos) const {
0098   vector<TransientTrack> ttVect;
0099   ttVect.reserve((*trkColl).size());
0100   for (unsigned int i = 0; i < (*trkColl).size(); i++) {
0101     TrackRef ref(trkColl, i);
0102     double time = trackTimes[ref];
0103     double timeReso = trackTimeResos[ref];
0104     timeReso = (timeReso > 1e-6 ? timeReso
0105                                 : defaultInvalidTrackTimeReso);  // make the error much larger than the BS time width
0106     if (edm::isNotFinite(time)) {
0107       time = 0.0;
0108       timeReso = defaultInvalidTrackTimeReso;
0109     }
0110     ttVect.push_back(TransientTrack(ref, time, timeReso, theField, theTrackingGeometry));
0111   }
0112   return ttVect;
0113 }
0114 
0115 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<reco::GsfTrackCollection>& trkColl,
0116                                                     const edm::ValueMap<float>& trackTimes,
0117                                                     const edm::ValueMap<float>& trackTimeResos) const {
0118   vector<TransientTrack> ttVect;
0119   ttVect.reserve((*trkColl).size());
0120   for (unsigned int i = 0; i < (*trkColl).size(); i++) {
0121     GsfTrackRef ref(trkColl, i);
0122     double time = trackTimes[ref];
0123     double timeReso = trackTimeResos[ref];
0124     timeReso = (timeReso > 1e-6 ? timeReso
0125                                 : defaultInvalidTrackTimeReso);  // make the error much larger than the BS time width
0126     if (edm::isNotFinite(time)) {
0127       time = 0.0;
0128       timeReso = defaultInvalidTrackTimeReso;
0129     }
0130     ttVect.push_back(TransientTrack(new GsfTransientTrack(ref, time, timeReso, theField, theTrackingGeometry)));
0131   }
0132   return ttVect;
0133 }
0134 
0135 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<edm::View<Track> >& trkColl,
0136                                                     const edm::ValueMap<float>& trackTimes,
0137                                                     const edm::ValueMap<float>& trackTimeResos) const {
0138   vector<TransientTrack> ttVect;
0139   ttVect.reserve((*trkColl).size());
0140   for (unsigned int i = 0; i < (*trkColl).size(); i++) {
0141     const Track* trk = &(*trkColl)[i];
0142     const GsfTrack* gsfTrack = dynamic_cast<const GsfTrack*>(trk);
0143     if (gsfTrack) {
0144       GsfTrackRef ref = RefToBase<Track>(trkColl, i).castTo<GsfTrackRef>();
0145       double time = trackTimes[ref];
0146       double timeReso = trackTimeResos[ref];
0147       timeReso = (timeReso > 1e-6 ? timeReso
0148                                   : defaultInvalidTrackTimeReso);  // make the error much larger than the BS time width
0149       if (edm::isNotFinite(time)) {
0150         time = 0.0;
0151         timeReso = defaultInvalidTrackTimeReso;
0152       }
0153       ttVect.push_back(TransientTrack(new GsfTransientTrack(
0154           RefToBase<Track>(trkColl, i).castTo<GsfTrackRef>(), time, timeReso, theField, theTrackingGeometry)));
0155     } else {  // gsf
0156       TrackRef ref = RefToBase<Track>(trkColl, i).castTo<TrackRef>();
0157       double time = trackTimes[ref];
0158       double timeReso = trackTimeResos[ref];
0159       timeReso = (timeReso > 1e-6 ? timeReso
0160                                   : defaultInvalidTrackTimeReso);  // make the error much larger than the BS time width
0161       if (edm::isNotFinite(time)) {
0162         time = 0.0;
0163         timeReso = defaultInvalidTrackTimeReso;
0164       }
0165       ttVect.push_back(TransientTrack(
0166           RefToBase<Track>(trkColl, i).castTo<TrackRef>(), time, timeReso, theField, theTrackingGeometry));
0167     }
0168   }
0169   return ttVect;
0170 }
0171 
0172 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<reco::TrackCollection>& trkColl,
0173                                                     const reco::BeamSpot& beamSpot) const {
0174   vector<TransientTrack> ttVect = build(trkColl);
0175   for (unsigned int i = 0; i < ttVect.size(); i++) {
0176     ttVect[i].setBeamSpot(beamSpot);
0177   }
0178   return ttVect;
0179 }
0180 
0181 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<reco::GsfTrackCollection>& trkColl,
0182                                                     const reco::BeamSpot& beamSpot) const {
0183   vector<TransientTrack> ttVect = build(trkColl);
0184   for (unsigned int i = 0; i < ttVect.size(); i++) {
0185     ttVect[i].setBeamSpot(beamSpot);
0186   }
0187   return ttVect;
0188 }
0189 
0190 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<edm::View<Track> >& trkColl,
0191                                                     const reco::BeamSpot& beamSpot) const {
0192   vector<TransientTrack> ttVect = build(trkColl);
0193   for (unsigned int i = 0; i < ttVect.size(); i++) {
0194     ttVect[i].setBeamSpot(beamSpot);
0195   }
0196   return ttVect;
0197 }
0198 
0199 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<reco::TrackCollection>& trkColl,
0200                                                     const reco::BeamSpot& beamSpot,
0201                                                     const edm::ValueMap<float>& trackTimes,
0202                                                     const edm::ValueMap<float>& trackTimeResos) const {
0203   vector<TransientTrack> ttVect = build(trkColl, trackTimes, trackTimeResos);
0204   for (unsigned int i = 0; i < ttVect.size(); i++) {
0205     ttVect[i].setBeamSpot(beamSpot);
0206   }
0207   return ttVect;
0208 }
0209 
0210 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<reco::GsfTrackCollection>& trkColl,
0211                                                     const reco::BeamSpot& beamSpot,
0212                                                     const edm::ValueMap<float>& trackTimes,
0213                                                     const edm::ValueMap<float>& trackTimeResos) const {
0214   vector<TransientTrack> ttVect = build(trkColl, trackTimes, trackTimeResos);
0215   for (unsigned int i = 0; i < ttVect.size(); i++) {
0216     ttVect[i].setBeamSpot(beamSpot);
0217   }
0218   return ttVect;
0219 }
0220 
0221 vector<TransientTrack> TransientTrackBuilder::build(const edm::Handle<edm::View<Track> >& trkColl,
0222                                                     const reco::BeamSpot& beamSpot,
0223                                                     const edm::ValueMap<float>& trackTimes,
0224                                                     const edm::ValueMap<float>& trackTimeResos) const {
0225   vector<TransientTrack> ttVect = build(trkColl, trackTimes, trackTimeResos);
0226   for (unsigned int i = 0; i < ttVect.size(); i++) {
0227     ttVect[i].setBeamSpot(beamSpot);
0228   }
0229   return ttVect;
0230 }
0231 
0232 TransientTrack TransientTrackBuilder::build(const FreeTrajectoryState& fts) const {
0233   return TransientTrack(new TransientTrackFromFTS(fts));
0234 }