CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoHI/HiTracking/plugins/HIPixelMedianVtxProducer.cc

Go to the documentation of this file.
00001 #include "HIPixelMedianVtxProducer.h"
00002 
00003 #include "FWCore/Framework/interface/MakerMacros.h"
00004 #include "FWCore/Framework/interface/ESHandle.h"
00005 #include "FWCore/Framework/interface/Frameworkfwd.h"
00006 #include "FWCore/Framework/interface/Event.h"
00007 
00008 #include "DataFormats/VertexReco/interface/Vertex.h"
00009 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00010 
00011 #include <fstream>
00012 #include <iostream>
00013 #include <vector>
00014 #include <algorithm>
00015 
00016 #include "TROOT.h"
00017 #include "TH1F.h"
00018 #include "TF1.h"
00019 
00020 /*****************************************************************************/
00021 HIPixelMedianVtxProducer::HIPixelMedianVtxProducer(const edm::ParameterSet& ps) : 
00022   theTrackCollection(ps.getParameter<edm::InputTag>("TrackCollection")),
00023   thePtMin(ps.getParameter<double>("PtMin")),
00024   thePeakFindThresh(ps.getParameter<unsigned int>("PeakFindThreshold")),
00025   thePeakFindMaxZ(ps.getParameter<double>("PeakFindMaxZ")),
00026   thePeakFindBinning(ps.getParameter<int>("PeakFindBinsPerCm")),
00027   theFitThreshold(ps.getParameter<int>("FitThreshold")),
00028   theFitMaxZ(ps.getParameter<double>("FitMaxZ")),
00029   theFitBinning(ps.getParameter<int>("FitBinsPerCm"))
00030 {
00031   produces<reco::VertexCollection>();
00032 }
00033 
00034 /*****************************************************************************/
00035 void HIPixelMedianVtxProducer::produce
00036 (edm::Event& ev, const edm::EventSetup& es)
00037 {
00038   // Get pixel tracks
00039   edm::Handle<reco::TrackCollection> trackCollection;
00040   ev.getByLabel(theTrackCollection, trackCollection);
00041   const reco::TrackCollection tracks_ = *(trackCollection.product());
00042   
00043   // Select tracks above minimum pt
00044   std::vector<const reco::Track *> tracks;
00045   for (unsigned int i=0; i<tracks_.size(); i++) {
00046     if (tracks_[i].pt() <  thePtMin && std::fabs(tracks_[i].vz()) < 100000.) continue;
00047     reco::TrackRef recTrack(trackCollection, i);
00048     tracks.push_back( &(*recTrack));
00049   }
00050         
00051   LogTrace("MinBiasTracking")
00052     << " [VertexProducer] selected tracks: "
00053     << tracks.size() << " (out of " << tracks_.size()
00054     << ") above pt = " << thePtMin; 
00055         
00056   // Output vertex collection
00057   std::auto_ptr<reco::VertexCollection> vertices(new reco::VertexCollection);
00058 
00059   // No tracks -> return empty collection
00060   if(tracks.size() == 0) {
00061       ev.put(vertices);
00062       return;
00063   }
00064 
00065   // Sort tracks according to vertex z position
00066   std::sort(tracks.begin(), tracks.end(), ComparePairs());
00067   
00068   // Calculate median vz
00069   float med;
00070   if(tracks.size() % 2 == 0)
00071     med = (tracks[tracks.size()/2-1]->vz() + tracks[tracks.size()/2]->vz())/2;
00072   else
00073     med =  tracks[tracks.size()/2  ]->vz();
00074   
00075   LogTrace("MinBiasTracking")
00076     << "  [vertex position] median    = " << med << " cm";
00077   
00078   // In high multiplicity events, fit around most probable position
00079   if(tracks.size() > thePeakFindThresh) { 
00080           
00081     // Find maximum bin
00082     TH1F hmax("hmax","hmax",thePeakFindBinning*2.0*thePeakFindMaxZ,-1.0*thePeakFindMaxZ,thePeakFindMaxZ);
00083     
00084     for(std::vector<const reco::Track *>::const_iterator
00085           track = tracks.begin(); track!= tracks.end(); track++)
00086       if(fabs((*track)->vz()) < thePeakFindMaxZ)
00087         hmax.Fill((*track)->vz());
00088     
00089     int maxBin = hmax.GetMaximumBin();
00090     
00091     LogTrace("MinBiasTracking")
00092       << "  [vertex position] most prob = "
00093       << hmax.GetBinCenter(maxBin) << " cm";
00094     
00095     // Find 3-bin weighted average
00096     float num=0.0, denom=0.0;
00097     for(int i=-1; i<=1; i++) {
00098       num += hmax.GetBinContent(maxBin+i)*hmax.GetBinCenter(maxBin+i);
00099       denom += hmax.GetBinContent(maxBin+i); 
00100     }
00101 
00102     if(denom==0) 
00103     {
00104       reco::Vertex::Error err;
00105       err(2,2) = 0.1 * 0.1;
00106       reco::Vertex ver(reco::Vertex::Point(0,0,99999.),
00107                        err, 0, 0, 1);
00108       vertices->push_back(ver);
00109       ev.put(vertices);
00110       return;
00111     }
00112 
00113     float nBinAvg = num/denom;
00114 
00115     // Center fit at 3-bin weighted average around max bin
00116     med = nBinAvg; 
00117     
00118     LogTrace("MinBiasTracking")
00119       << "  [vertex position] 3-bin weighted average = "
00120       << nBinAvg << " cm";
00121   }
00122   
00123   // Bin vz-values around most probable value (or median) for fitting
00124   TH1F histo("histo","histo", theFitBinning*2.0*theFitMaxZ,-1.0*theFitMaxZ,theFitMaxZ);
00125   histo.Sumw2();
00126 
00127   for(std::vector<const reco::Track *>::const_iterator
00128         track = tracks.begin(); track!= tracks.end(); track++)
00129     if(fabs((*track)->vz() - med) < theFitMaxZ)
00130       histo.Fill((*track)->vz() - med);
00131   
00132   LogTrace("MinBiasTracking")
00133     << "  [vertex position] most prob for fit = "
00134     << med + histo.GetBinCenter(histo.GetMaximumBin()) << " cm";
00135 
00136   // If there are very few entries, don't do the fit
00137   if(histo.GetEntries() <= theFitThreshold) {
00138 
00139     LogTrace("MinBiasTracking")
00140       << "  [vertex position] Fewer than" << theFitThreshold
00141       <<  " entries in fit histogram. Returning median.";
00142     
00143     reco::Vertex::Error err;
00144     err(2,2) = 0.1 * 0.1;
00145     reco::Vertex ver(reco::Vertex::Point(0,0,med),
00146                      err, 0, 1, 1);
00147     vertices->push_back(ver);
00148     ev.put(vertices);
00149     return;
00150   }
00151 
00152   // Otherwise, there are enough entries to refine the estimate with a fit
00153   TF1 f1("f1","[0]*exp(-0.5 * ((x-[1])/[2])^2) + [3]");
00154   f1.SetParameters(10.,0.,0.02,0.002*tracks.size());
00155   f1.SetParLimits(1,-0.1,0.1);
00156   f1.SetParLimits(2,0.001,0.05);
00157   f1.SetParLimits(3,0.0,0.005*tracks.size());
00158     
00159   histo.Fit("f1","QN");
00160     
00161   LogTrace("MinBiasTracking")
00162     << "  [vertex position] fitted    = "
00163     << med + f1.GetParameter(1) << " +- " << f1.GetParError(1) << " cm";
00164   
00165   reco::Vertex::Error err;
00166   err(2,2) = f1.GetParError(1) * f1.GetParError(1); 
00167   reco::Vertex ver(reco::Vertex::Point(0,0,med + f1.GetParameter(1)),
00168                    err, 0, 1, 1);
00169   vertices->push_back(ver);
00170   
00171   ev.put(vertices);
00172   return;
00173 
00174 }
00175