#include <PixelVertexProducerMedian.h>
Public Member Functions | |
PixelVertexProducerMedian (const edm::ParameterSet &ps) | |
virtual void | produce (edm::Event &ev, const edm::EventSetup &es) |
~PixelVertexProducerMedian () | |
Private Attributes | |
edm::ParameterSet | theConfig |
double | thePtMin |
Definition at line 9 of file PixelVertexProducerMedian.h.
PixelVertexProducerMedian::PixelVertexProducerMedian | ( | const edm::ParameterSet & | ps | ) | [explicit] |
Definition at line 35 of file PixelVertexProducerMedian.cc.
: theConfig(ps) { produces<reco::VertexCollection>(); }
PixelVertexProducerMedian::~PixelVertexProducerMedian | ( | ) |
Definition at line 42 of file PixelVertexProducerMedian.cc.
{ }
void PixelVertexProducerMedian::produce | ( | edm::Event & | ev, |
const edm::EventSetup & | es | ||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 48 of file PixelVertexProducerMedian.cc.
References validate-o2o-wbm::f1, edm::Event::getByLabel(), timingPdfMaker::histo, i, LogTrace, edm::Handle< T >::product(), edm::Event::put(), python::multivaluedict::sort(), testEve_cfg::tracks, and tracks_.
{ // Get pixel tracks edm::Handle<reco::TrackCollection> trackCollection; std::string trackCollectionName = theConfig.getParameter<std::string>("TrackCollection"); ev.getByLabel(trackCollectionName, trackCollection); const reco::TrackCollection tracks_ = *(trackCollection.product()); thePtMin = theConfig.getParameter<double>("PtMin"); // Select tracks std::vector<const reco::Track *> tracks; for (unsigned int i=0; i<tracks_.size(); i++) { if (tracks_[i].pt() > thePtMin) { reco::TrackRef recTrack(trackCollection, i); tracks.push_back( &(*recTrack)); } } LogTrace("MinBiasTracking") << " [VertexProducer] selected tracks: " << tracks.size() << " (out of " << tracks_.size() << ")"; std::auto_ptr<reco::VertexCollection> vertices(new reco::VertexCollection); if(tracks.size() > 0) { // Sort along vertex z position std::sort(tracks.begin(), tracks.end(), ComparePairs()); // Median float med; if(tracks.size() % 2 == 0) med = (tracks[tracks.size()/2-1]->vz() + tracks[tracks.size()/2]->vz())/2; else med = tracks[tracks.size()/2 ]->vz(); LogTrace("MinBiasTracking") << " [vertex position] median = " << med << " cm"; if(tracks.size() > 10) { // Binning around med, halfWidth int nBin = 100; float halfWidth = 0.1; // cm // Most probable TH1F histo("histo","histo", nBin, -halfWidth,halfWidth); for(std::vector<const reco::Track *>::const_iterator track = tracks.begin(); track!= tracks.end(); track++) if(fabs((*track)->vz() - med) < halfWidth) histo.Fill((*track)->vz() - med); LogTrace("MinBiasTracking") << " [vertex position] most prob = " << med + histo.GetBinCenter(histo.GetMaximumBin()) << " cm"; // Fit above max/2 histo.Sumw2(); TF1 f1("f1","[0]*exp(-0.5 * ((x-[1])/[2])^2) + [3]"); f1.SetParameters(10.,0.,0.01, 1.); histo.Fit("f1","QN"); LogTrace("MinBiasTracking") << " [vertex position] fitted = " << med + f1.GetParameter(1) << " +- " << f1.GetParError(1) << " cm"; // Store reco::Vertex::Error err; err(2,2) = f1.GetParError(1) * f1.GetParError(1); reco::Vertex ver(reco::Vertex::Point(0,0,med + f1.GetParameter(1)), err, 0, 1, 1); vertices->push_back(ver); } else { // Store reco::Vertex::Error err; err(2,2) = 0.1 * 0.1; reco::Vertex ver(reco::Vertex::Point(0,0,med), err, 0, 1, 1); vertices->push_back(ver); } } ev.put(vertices); }
Definition at line 17 of file PixelVertexProducerMedian.h.
double PixelVertexProducerMedian::thePtMin [private] |
Definition at line 18 of file PixelVertexProducerMedian.h.