CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoTracker/TkSeedGenerator/plugins/SeedGeneratorFromProtoTracksEDProducer.cc

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/VertexReco/interface/Vertex.h"
00008 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00009 
00010 #include "RecoTracker/TkSeedGenerator/interface/SeedFromProtoTrack.h"
00011 #include "RecoTracker/TkSeedingLayers/interface/SeedingHitSet.h"
00012 #include "SeedFromConsecutiveHitsCreator.h"
00013 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00014 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00015 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
00016 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00017 
00018 #include "FWCore/Framework/interface/EventSetup.h"
00019 #include "FWCore/Framework/interface/ESHandle.h"
00020 
00021 
00022 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00023 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025 #include <vector>
00026 #include <cassert>
00027 using namespace edm;
00028 using namespace reco;
00029 
00030 template <class T> T sqr( T t) {return t*t;}
00031 typedef TransientTrackingRecHit::ConstRecHitPointer Hit;
00032 
00033 struct HitLessByRadius { bool operator() (const Hit& h1, const Hit & h2) { return h1->globalPosition().perp2() < h2->globalPosition().perp2(); } };
00034 
00035 SeedGeneratorFromProtoTracksEDProducer::SeedGeneratorFromProtoTracksEDProducer(const ParameterSet& cfg)
00036   : theConfig(cfg), 
00037     theInputCollectionTag(cfg.getParameter<InputTag>("InputCollection")),
00038     theInputVertexCollectionTag(cfg.getParameter<InputTag>("InputVertexCollection")),
00039     originHalfLength(cfg.getParameter<double>("originHalfLength")),
00040     originRadius(cfg.getParameter<double>("originRadius")),
00041     useProtoTrackKinematics(cfg.getParameter<bool>("useProtoTrackKinematics")),
00042     useEventsWithNoVertex(cfg.getParameter<bool>("useEventsWithNoVertex")),
00043     builderName(cfg.getParameter<std::string>("TTRHBuilder"))
00044 
00045 {
00046   produces<TrajectorySeedCollection>();
00047 }
00048 
00049 void SeedGeneratorFromProtoTracksEDProducer::produce(edm::Event& ev, const edm::EventSetup& es)
00050 {
00051   std::auto_ptr<TrajectorySeedCollection> result(new TrajectorySeedCollection());
00052   Handle<reco::TrackCollection> trks;
00053   ev.getByLabel(theInputCollectionTag, trks);
00054 
00055   const TrackCollection &protos = *(trks.product());
00056   
00057   edm::Handle<reco::VertexCollection> vertices;
00058   bool foundVertices = ev.getByLabel(theInputVertexCollectionTag, vertices);
00059   //const reco::VertexCollection & vertices = *(h_vertices.product());
00060 
00061   for (TrackCollection::const_iterator it=protos.begin(); it!= protos.end(); ++it) {
00062     const Track & proto = (*it);
00063     GlobalPoint vtx(proto.vertex().x(), proto.vertex().y(), proto.vertex().z());
00064 
00065     // check the compatibility with a primary vertex
00066     bool keepTrack = false;
00067     if ( !foundVertices ) { 
00068           if (useEventsWithNoVertex) keepTrack = true;
00069     } else { 
00070       for (reco::VertexCollection::const_iterator iv=vertices->begin(); iv!= vertices->end(); ++iv) {
00071         GlobalPoint aPV(iv->position().x(),iv->position().y(),iv->position().z());
00072         double distR2 = sqr(vtx.x()-aPV.x()) +sqr(vtx.y()-aPV.y());
00073         double distZ = fabs(vtx.z()-aPV.z());
00074         if ( distR2 < sqr(originRadius) && distZ < originHalfLength ) { 
00075           keepTrack = true;
00076           break;
00077         }
00078       }
00079     }
00080     if (!keepTrack) continue;
00081 
00082     if ( useProtoTrackKinematics ) {
00083       SeedFromProtoTrack seedFromProtoTrack( proto, es);
00084       if (seedFromProtoTrack.isValid()) (*result).push_back( seedFromProtoTrack.trajectorySeed() );
00085     } else {
00086       edm::ESHandle<TransientTrackingRecHitBuilder> ttrhbESH;
00087       es.get<TransientRecHitRecord>().get(builderName,ttrhbESH);
00088       std::vector<Hit> hits;
00089       for (unsigned int iHit = 0, nHits = proto.recHitsSize(); iHit < nHits; ++iHit) {
00090         TrackingRecHitRef refHit = proto.recHit(iHit);
00091         if(refHit->isValid()) hits.push_back(ttrhbESH->build(  &(*refHit) ));
00092       }
00093       sort(hits.begin(), hits.end(), HitLessByRadius());
00094       assert(hits.size()<4);
00095       if (hits.size() > 1) {
00096         double mom_perp = sqrt(proto.momentum().x()*proto.momentum().x()+proto.momentum().y()*proto.momentum().y());
00097         GlobalTrackingRegion region(mom_perp, vtx, 0.2, 0.2);
00098         SeedFromConsecutiveHitsCreator().trajectorySeed(*result, SeedingHitSet(hits[0], hits[1], hits.size() >2 ? hits[2] : SeedingHitSet::nullPtr() ), region, es, 0);
00099       }
00100     }
00101   } 
00102 
00103   ev.put(result);
00104 }
00105