Go to the documentation of this file.00001 #include "SeedGeneratorFromProtoTracksEDProducer.h"
00002
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "DataFormats/Common/interface/Handle.h"
00005
00006 #include "DataFormats/TrackReco/interface/Track.h"
00007 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00008
00009 #include "RecoTracker/TkSeedGenerator/interface/SeedFromProtoTrack.h"
00010 #include "RecoTracker/TkSeedingLayers/interface/SeedingHitSet.h"
00011 #include "RecoTracker/TkSeedGenerator/interface/SeedFromConsecutiveHitsCreator.h"
00012 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00013 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00014 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
00015 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00016
00017 #include "FWCore/Framework/interface/EventSetup.h"
00018 #include "FWCore/Framework/interface/ESHandle.h"
00019
00020
00021 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 #include <vector>
00024
00025 using namespace edm;
00026 using namespace reco;
00027
00028 template <class T> T sqr( T t) {return t*t;}
00029 typedef TransientTrackingRecHit::ConstRecHitPointer Hit;
00030
00031 struct HitLessByRadius { bool operator() (const Hit& h1, const Hit & h2) { return h1->globalPosition().perp2() < h2->globalPosition().perp2(); } };
00032
00033 SeedGeneratorFromProtoTracksEDProducer::SeedGeneratorFromProtoTracksEDProducer(const ParameterSet& cfg)
00034 : theConfig(cfg), theInputCollectionTag(cfg.getParameter<InputTag>("InputCollection"))
00035 {
00036 produces<TrajectorySeedCollection>();
00037 }
00038
00039 void SeedGeneratorFromProtoTracksEDProducer::produce(edm::Event& ev, const edm::EventSetup& es)
00040 {
00041 std::auto_ptr<TrajectorySeedCollection> result(new TrajectorySeedCollection());
00042 Handle<reco::TrackCollection> trks;
00043 ev.getByLabel(theInputCollectionTag, trks);
00044
00045 const TrackCollection &protos = *(trks.product());
00046
00047 for (TrackCollection::const_iterator it=protos.begin(); it!= protos.end(); ++it) {
00048 const Track & proto = (*it);
00049
00050 if (theConfig.getParameter<bool>("useProtoTrackKinematics")) {
00051 SeedFromProtoTrack seedFromProtoTrack( proto, es);
00052 if (seedFromProtoTrack.isValid()) (*result).push_back( seedFromProtoTrack.trajectorySeed() );
00053 } else {
00054 edm::ESHandle<TransientTrackingRecHitBuilder> ttrhbESH;
00055 std::string builderName = theConfig.getParameter<std::string>("TTRHBuilder");
00056 es.get<TransientRecHitRecord>().get(builderName,ttrhbESH);
00057 std::vector<Hit> hits;
00058 for (unsigned int iHit = 0, nHits = proto.recHitsSize(); iHit < nHits; ++iHit) {
00059 TrackingRecHitRef refHit = proto.recHit(iHit);
00060 if(refHit->isValid()) hits.push_back(ttrhbESH->build( &(*refHit) ));
00061 sort(hits.begin(), hits.end(), HitLessByRadius());
00062 }
00063 if (hits.size() >= 2) {
00064 GlobalPoint vtx(proto.vertex().x(), proto.vertex().y(), proto.vertex().z());
00065 double mom_perp = sqrt(proto.momentum().x()*proto.momentum().x()+proto.momentum().y()*proto.momentum().y());
00066 GlobalTrackingRegion region(mom_perp, vtx, 0.2, 0.2);
00067 SeedFromConsecutiveHitsCreator().trajectorySeed(*result, SeedingHitSet(hits), region, es);
00068 }
00069 }
00070 }
00071
00072 ev.put(result);
00073 }
00074