![]() |
![]() |
#include <RecoPixelVertexing/PixelLowPtUtilities/plugins/PixelTrackProducerWithZPos.h>
Public Member Functions | |
PixelTrackProducerWithZPos (const edm::ParameterSet &conf) | |
virtual void | produce (edm::Event &ev, const edm::EventSetup &es) |
~PixelTrackProducerWithZPos () | |
Private Member Functions | |
void | beginJob (const edm::EventSetup &es) |
std::pair< float, float > | refitWithVertex (const reco::Track &recTrack, const reco::VertexCollection *vertices) |
void | store (edm::Event &ev, const pixeltrackfitting::TracksWithRecHits &selectedTracks) |
Private Attributes | |
edm::ParameterSet | ps |
PixelTrackCleaner * | theCleaner |
const TrackHitsFilter * | theFilter |
const PixelFitter * | theFitter |
OrderedHitsGenerator * | theGenerator |
const TrackHitsFilter * | theHitsFilter |
double | theOriginRadius |
double | thePtMin |
TrackingRegionProducer * | theRegionProducer |
const TransientTrackBuilder * | theTTBuilder |
bool | theUseFoundVertices |
Definition at line 26 of file PixelTrackProducerWithZPos.h.
PixelTrackProducerWithZPos::PixelTrackProducerWithZPos | ( | const edm::ParameterSet & | conf | ) | [explicit] |
Definition at line 75 of file PixelTrackProducerWithZPos.cc.
00076 : ps(conf), theFitter(0), theFilter(0), theHitsFilter(0), theCleaner(0), theGenerator(0), theRegionProducer(0) 00077 { 00078 edm::LogInfo("PixelTrackProducerWithZPos")<<" construction..."; 00079 produces<reco::TrackCollection>(); 00080 produces<TrackingRecHitCollection>(); 00081 produces<reco::TrackExtraCollection>(); 00082 }
PixelTrackProducerWithZPos::~PixelTrackProducerWithZPos | ( | ) |
Definition at line 86 of file PixelTrackProducerWithZPos.cc.
References theCleaner, theFilter, theFitter, theGenerator, theHitsFilter, and theRegionProducer.
00087 { 00088 delete theFilter; 00089 delete theHitsFilter; 00090 delete theFitter; 00091 delete theCleaner; 00092 delete theGenerator; 00093 delete theRegionProducer; 00094 }
void PixelTrackProducerWithZPos::beginJob | ( | const edm::EventSetup & | es | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 97 of file PixelTrackProducerWithZPos.cc.
References edm::EventSetup::get(), DBSPlugin::get(), edm::ParameterSet::getParameter(), edm::ESHandle< T >::product(), ps, theCleaner, theFilter, theFitter, theGenerator, theHitsFilter, theOriginRadius, thePtMin, theRegionProducer, theTTBuilder, and theUseFoundVertices.
00098 { 00099 // Region 00100 ParameterSet regfactoryPSet = ps.getParameter<ParameterSet>("RegionFactoryPSet"); 00101 std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName"); 00102 theRegionProducer = TrackingRegionProducerFactory::get()->create(regfactoryName,regfactoryPSet); 00103 00104 // OrderesHits 00105 ParameterSet orderedPSet = ps.getParameter<ParameterSet>("OrderedHitsFactoryPSet"); 00106 std::string orderedName = orderedPSet.getParameter<std::string>("ComponentName"); 00107 theGenerator = OrderedHitsGeneratorFactory::get()->create( orderedName, orderedPSet); 00108 00109 // Fitter 00110 ParameterSet fitterPSet = ps.getParameter<ParameterSet>("FitterPSet"); 00111 std::string fitterName = fitterPSet.getParameter<std::string>("ComponentName"); 00112 theFitter = PixelFitterFactory::get()->create( fitterName, fitterPSet); 00113 00114 // Filter 00115 ParameterSet filterPSet = ps.getParameter<ParameterSet>("FilterPSet"); 00116 std::string filterName = filterPSet.getParameter<std::string>("ComponentName"); 00117 theFilter = TrackHitsFilterFactory::get()->create( filterName, filterPSet, es); 00118 00119 // Cleaner 00120 theHitsFilter = TrackHitsFilterFactory::get()->create("ValidHitPairFilter", filterPSet, es); 00121 00122 ParameterSet cleanerPSet = ps.getParameter<ParameterSet>("CleanerPSet"); 00123 std::string cleanerName = cleanerPSet.getParameter<std::string>("ComponentName"); 00124 theCleaner = PixelTrackCleanerFactory::get()->create( cleanerName, cleanerPSet); 00125 00126 // Get transient track builder 00127 edm::ParameterSet regionPSet = regfactoryPSet.getParameter<edm::ParameterSet>("RegionPSet"); 00128 theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices"); 00129 // theUseChi2Cut = regionPSet.getParameter<bool>("useChi2Cut"); 00130 thePtMin = regionPSet.getParameter<double>("ptMin"); 00131 theOriginRadius = regionPSet.getParameter<double>("originRadius"); 00132 00133 if(theUseFoundVertices) 00134 { 00135 edm::ESHandle<TransientTrackBuilder> builder; 00136 es.get<TransientTrackRecord>().get("TransientTrackBuilder", builder); 00137 theTTBuilder = builder.product(); 00138 } 00139 }
void PixelTrackProducerWithZPos::produce | ( | edm::Event & | ev, | |
const edm::EventSetup & | es | |||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 188 of file PixelTrackProducerWithZPos.cc.
References PixelTrackCleaner::cleanTracks(), edm::Event::getByType(), edm::ParameterSet::getParameter(), LogTrace, edm::Handle< T >::product(), ps, TrackingRegionProducer::regions(), PixelFitter::run(), OrderedHitsGenerator::run(), SeedingHitSet::size(), OrderedSeedingHits::size(), store(), theCleaner, theFilter, theFitter, theGenerator, theHitsFilter, theRegionProducer, theUseFoundVertices, track, and tracks.
00189 { 00190 LogTrace("MinBiasTracking") 00191 << "\033[22;31m[" 00192 << ps.getParameter<string>("passLabel") 00193 << "]\033[22;0m"; 00194 00195 TracksWithRecHits tracks; 00196 00197 // Get vertices 00198 const reco::VertexCollection* vertices = 0; 00199 00200 if(theUseFoundVertices) 00201 { 00202 edm::Handle<reco::VertexCollection> vertexCollection; 00203 ev.getByType(vertexCollection); 00204 vertices = vertexCollection.product(); 00205 } 00206 00207 typedef std::vector<TrackingRegion* > Regions; 00208 typedef Regions::const_iterator IR; 00209 Regions regions = theRegionProducer->regions(ev,es); 00210 00211 for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) 00212 { 00213 const TrackingRegion & region = **ir; 00214 00215 const OrderedSeedingHits & triplets = theGenerator->run(region,ev,es); 00216 unsigned int nTriplets = triplets.size(); 00217 00218 LogTrace("MinBiasTracking") 00219 << " [TrackProducer] number of triplets : " 00220 << triplets.size(); 00221 00222 // producing tracks 00223 for(unsigned int iTriplet = 0; iTriplet < nTriplets; ++iTriplet) 00224 { 00225 const SeedingHitSet & triplet = triplets[iTriplet]; 00226 00227 std::vector<const TrackingRecHit *> hits; 00228 for (unsigned int iHit = 0, nHits = triplet.size(); iHit < nHits; ++iHit) 00229 hits.push_back( triplet[iHit] ); 00230 00231 // Fitter 00232 reco::Track* track = theFitter->run(es, hits, region); 00233 00234 if(hits.size() == 2) 00235 if ( ! (*theHitsFilter)(track, hits) ) 00236 { 00237 delete track; 00238 continue; 00239 } 00240 00241 // Filter 00242 if ( ! (*theFilter)(track,hits) ) 00243 { 00244 LogTrace("MinBiasTracking") << " [TrackProducer] track did not pass cluster shape filter"; 00245 delete track; 00246 continue; 00247 } 00248 00249 // add tracks 00250 tracks.push_back(TrackWithRecHits(track, hits)); 00251 } 00252 } 00253 00254 // Cleaner 00255 if(theCleaner) tracks = theCleaner->cleanTracks(tracks); 00256 00257 // store tracks 00258 store(ev, tracks); 00259 00260 // clean memory 00261 for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) 00262 delete (*ir); 00263 }
pair< float, float > PixelTrackProducerWithZPos::refitWithVertex | ( | const reco::Track & | recTrack, | |
const reco::VertexCollection * | vertices | |||
) | [private] |
Definition at line 143 of file PixelTrackProducerWithZPos.cc.
References TransientTrackBuilder::build(), reco::Vertex::covariance(), e, reco::Vertex::position(), reco::TrackBase::pt(), HLT_VtxMuL3::result, theTTBuilder, and reco::TrackBase::vertex().
00145 { 00146 TransientTrack theTransientTrack = theTTBuilder->build(recTrack); 00147 00148 // If there are vertices found 00149 if(vertices->size() > 0) 00150 { 00151 float dzmin = -1.; 00152 const reco::Vertex * closestVertex = 0; 00153 00154 // Look for the closest vertex in z 00155 for(reco::VertexCollection::const_iterator 00156 vertex = vertices->begin(); vertex!= vertices->end(); vertex++) 00157 { 00158 float dz = fabs(recTrack.vertex().z() - vertex->position().z()); 00159 if(vertex == vertices->begin() || dz < dzmin) 00160 { dzmin = dz ; closestVertex = &(*vertex); } 00161 } 00162 00163 00164 // Get vertex position and error matrix 00165 GlobalPoint vertexPosition(closestVertex->position().x(), 00166 closestVertex->position().y(), 00167 closestVertex->position().z()); 00168 00169 float beamSize = 15e-4; // 15 um 00170 GlobalError vertexError(beamSize*beamSize, 0, 00171 beamSize*beamSize, 0, 00172 0,closestVertex->covariance(2,2)); 00173 00174 // Refit track with vertex constraint 00175 SingleTrackVertexConstraint stvc; 00176 pair<TransientTrack, float> result = 00177 stvc.constrain(theTransientTrack, vertexPosition, vertexError); 00178 00179 return pair<float,float>(result.first.impactPointTSCP().pt(), 00180 result.second); 00181 } 00182 else 00183 return pair<float,float>(recTrack.pt(), -9999); 00184 }
void PixelTrackProducerWithZPos::store | ( | edm::Event & | ev, | |
const pixeltrackfitting::TracksWithRecHits & | selectedTracks | |||
) | [private] |
Definition at line 267 of file PixelTrackProducerWithZPos.cc.
References reco::TrackExtraBase::add(), clone(), i, k, LogDebug, edm::Event::put(), reco::TrackBase::setHitPattern(), track, and tracks.
Referenced by produce().
00268 { 00269 std::auto_ptr<reco::TrackCollection> tracks(new reco::TrackCollection); 00270 std::auto_ptr<TrackingRecHitCollection> recHits(new TrackingRecHitCollection); 00271 std::auto_ptr<reco::TrackExtraCollection> trackExtras(new reco::TrackExtraCollection); 00272 typedef std::vector<const TrackingRecHit *> RecHits; 00273 00274 int cc = 0, nTracks = cleanedTracks.size(); 00275 00276 for (int i = 0; i < nTracks; i++) 00277 { 00278 reco::Track* track = cleanedTracks.at(i).first; 00279 const RecHits & hits = cleanedTracks.at(i).second; 00280 00281 for (unsigned int k = 0; k < hits.size(); k++) 00282 { 00283 TrackingRecHit *hit = (hits.at(k))->clone(); 00284 track->setHitPattern(*hit, k); 00285 recHits->push_back(hit); 00286 } 00287 tracks->push_back(*track); 00288 delete track; 00289 00290 } 00291 00292 LogDebug("TrackProducer") << "put the collection of TrackingRecHit in the event" << "\n"; 00293 edm::OrphanHandle <TrackingRecHitCollection> ohRH = ev.put( recHits ); 00294 00295 00296 for (int k = 0; k < nTracks; k++) 00297 { 00298 reco::TrackExtra* theTrackExtra = new reco::TrackExtra(); 00299 00300 //fill the TrackExtra with TrackingRecHitRef 00301 unsigned int nHits = tracks->at(k).numberOfValidHits(); 00302 for(unsigned int i = 0; i < nHits; ++i) { 00303 theTrackExtra->add(TrackingRecHitRef(ohRH,cc)); 00304 cc++; 00305 } 00306 00307 trackExtras->push_back(*theTrackExtra); 00308 delete theTrackExtra; 00309 } 00310 00311 LogDebug("TrackProducer") << "put the collection of TrackExtra in the event" << "\n"; 00312 edm::OrphanHandle<reco::TrackExtraCollection> ohTE = ev.put(trackExtras); 00313 00314 for (int k = 0; k < nTracks; k++) 00315 { 00316 const reco::TrackExtraRef theTrackExtraRef(ohTE,k); 00317 (tracks->at(k)).setExtra(theTrackExtraRef); 00318 } 00319 00320 ev.put(tracks); 00321 }
Definition at line 41 of file PixelTrackProducerWithZPos.h.
Referenced by beginJob(), and produce().
Definition at line 46 of file PixelTrackProducerWithZPos.h.
Referenced by beginJob(), produce(), and ~PixelTrackProducerWithZPos().
const TrackHitsFilter* PixelTrackProducerWithZPos::theFilter [private] |
Definition at line 44 of file PixelTrackProducerWithZPos.h.
Referenced by beginJob(), produce(), and ~PixelTrackProducerWithZPos().
const PixelFitter* PixelTrackProducerWithZPos::theFitter [private] |
Definition at line 43 of file PixelTrackProducerWithZPos.h.
Referenced by beginJob(), produce(), and ~PixelTrackProducerWithZPos().
Definition at line 47 of file PixelTrackProducerWithZPos.h.
Referenced by beginJob(), produce(), and ~PixelTrackProducerWithZPos().
const TrackHitsFilter* PixelTrackProducerWithZPos::theHitsFilter [private] |
Definition at line 45 of file PixelTrackProducerWithZPos.h.
Referenced by beginJob(), produce(), and ~PixelTrackProducerWithZPos().
double PixelTrackProducerWithZPos::theOriginRadius [private] |
double PixelTrackProducerWithZPos::thePtMin [private] |
Definition at line 48 of file PixelTrackProducerWithZPos.h.
Referenced by beginJob(), produce(), and ~PixelTrackProducerWithZPos().
const TransientTrackBuilder* PixelTrackProducerWithZPos::theTTBuilder [private] |
Definition at line 50 of file PixelTrackProducerWithZPos.h.
Referenced by beginJob(), and refitWithVertex().
Definition at line 51 of file PixelTrackProducerWithZPos.h.
Referenced by beginJob(), and produce().