CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HIPixelMedianVtxProducer.cc
Go to the documentation of this file.
2 
7 
10 
11 #include <fstream>
12 #include <iostream>
13 #include <vector>
14 #include <algorithm>
15 
16 #include "TROOT.h"
17 #include "TH1F.h"
18 #include "TF1.h"
19 
20 /*****************************************************************************/
22  theTrackCollection(consumes<reco::TrackCollection>(ps.getParameter<edm::InputTag>("TrackCollection"))),
23  thePtMin(ps.getParameter<double>("PtMin")),
24  thePeakFindThresh(ps.getParameter<unsigned int>("PeakFindThreshold")),
25  thePeakFindMaxZ(ps.getParameter<double>("PeakFindMaxZ")),
26  thePeakFindBinning(ps.getParameter<int>("PeakFindBinsPerCm")),
27  theFitThreshold(ps.getParameter<int>("FitThreshold")),
28  theFitMaxZ(ps.getParameter<double>("FitMaxZ")),
29  theFitBinning(ps.getParameter<int>("FitBinsPerCm"))
30 {
31  produces<reco::VertexCollection>();
32 }
33 
34 /*****************************************************************************/
37 {
38  // Get pixel tracks
40  ev.getByToken(theTrackCollection, trackCollection);
41  const reco::TrackCollection tracks_ = *(trackCollection.product());
42 
43  // Select tracks above minimum pt
44  std::vector<const reco::Track *> tracks;
45  for (unsigned int i=0; i<tracks_.size(); i++) {
46  if (tracks_[i].pt() < thePtMin && std::fabs(tracks_[i].vz()) < 100000.) continue;
47  reco::TrackRef recTrack(trackCollection, i);
48  tracks.push_back( &(*recTrack));
49  }
50 
51  LogTrace("MinBiasTracking")
52  << " [VertexProducer] selected tracks: "
53  << tracks.size() << " (out of " << tracks_.size()
54  << ") above pt = " << thePtMin;
55 
56  // Output vertex collection
57  std::auto_ptr<reco::VertexCollection> vertices(new reco::VertexCollection);
58 
59  // No tracks -> return empty collection
60  if(tracks.size() == 0) {
61  ev.put(vertices);
62  return;
63  }
64 
65  // Sort tracks according to vertex z position
66  std::sort(tracks.begin(), tracks.end(), ComparePairs());
67 
68  // Calculate median vz
69  float med;
70  if(tracks.size() % 2 == 0)
71  med = (tracks[tracks.size()/2-1]->vz() + tracks[tracks.size()/2]->vz())/2;
72  else
73  med = tracks[tracks.size()/2 ]->vz();
74 
75  LogTrace("MinBiasTracking")
76  << " [vertex position] median = " << med << " cm";
77 
78  // In high multiplicity events, fit around most probable position
79  if(tracks.size() > thePeakFindThresh) {
80 
81  // Find maximum bin
82  TH1F hmax("hmax","hmax",thePeakFindBinning*2.0*thePeakFindMaxZ,-1.0*thePeakFindMaxZ,thePeakFindMaxZ);
83 
84  for(std::vector<const reco::Track *>::const_iterator
85  track = tracks.begin(); track!= tracks.end(); track++)
86  if(fabs((*track)->vz()) < thePeakFindMaxZ)
87  hmax.Fill((*track)->vz());
88 
89  int maxBin = hmax.GetMaximumBin();
90 
91  LogTrace("MinBiasTracking")
92  << " [vertex position] most prob = "
93  << hmax.GetBinCenter(maxBin) << " cm";
94 
95  // Find 3-bin weighted average
96  float num=0.0, denom=0.0;
97  for(int i=-1; i<=1; i++) {
98  num += hmax.GetBinContent(maxBin+i)*hmax.GetBinCenter(maxBin+i);
99  denom += hmax.GetBinContent(maxBin+i);
100  }
101 
102  if(denom==0)
103  {
105  err(2,2) = 0.1 * 0.1;
106  reco::Vertex ver(reco::Vertex::Point(0,0,99999.),
107  err, 0, 0, 1);
108  vertices->push_back(ver);
109  ev.put(vertices);
110  return;
111  }
112 
113  float nBinAvg = num/denom;
114 
115  // Center fit at 3-bin weighted average around max bin
116  med = nBinAvg;
117 
118  LogTrace("MinBiasTracking")
119  << " [vertex position] 3-bin weighted average = "
120  << nBinAvg << " cm";
121  }
122 
123  // Bin vz-values around most probable value (or median) for fitting
124  TH1F histo("histo","histo", theFitBinning*2.0*theFitMaxZ,-1.0*theFitMaxZ,theFitMaxZ);
125  histo.Sumw2();
126 
127  for(std::vector<const reco::Track *>::const_iterator
128  track = tracks.begin(); track!= tracks.end(); track++)
129  if(fabs((*track)->vz() - med) < theFitMaxZ)
130  histo.Fill((*track)->vz() - med);
131 
132  LogTrace("MinBiasTracking")
133  << " [vertex position] most prob for fit = "
134  << med + histo.GetBinCenter(histo.GetMaximumBin()) << " cm";
135 
136  // If there are very few entries, don't do the fit
137  if(histo.GetEntries() <= theFitThreshold) {
138 
139  LogTrace("MinBiasTracking")
140  << " [vertex position] Fewer than" << theFitThreshold
141  << " entries in fit histogram. Returning median.";
142 
144  err(2,2) = 0.1 * 0.1;
145  reco::Vertex ver(reco::Vertex::Point(0,0,med),
146  err, 0, 1, 1);
147  vertices->push_back(ver);
148  ev.put(vertices);
149  return;
150  }
151 
152  // Otherwise, there are enough entries to refine the estimate with a fit
153  TF1 f1("f1","[0]*exp(-0.5 * ((x-[1])/[2])^2) + [3]");
154  f1.SetParameters(10.,0.,0.02,0.002*tracks.size());
155  f1.SetParLimits(1,-0.1,0.1);
156  f1.SetParLimits(2,0.001,0.05);
157  f1.SetParLimits(3,0.0,0.005*tracks.size());
158 
159  histo.Fit(&f1,"QN");
160 
161  LogTrace("MinBiasTracking")
162  << " [vertex position] fitted = "
163  << med + f1.GetParameter(1) << " +- " << f1.GetParError(1) << " cm";
164 
166  err(2,2) = f1.GetParError(1) * f1.GetParError(1);
167  reco::Vertex ver(reco::Vertex::Point(0,0,med + f1.GetParameter(1)),
168  err, 0, 1, 1);
169  vertices->push_back(ver);
170 
171  ev.put(vertices);
172  return;
173 
174 }
175 
int i
Definition: DBlmapReader.cc:9
const std::vector< reco::PFCandidatePtr > & tracks_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
virtual void produce(edm::Event &ev, const edm::EventSetup &es)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
#define LogTrace(id)
T const * product() const
Definition: Handle.h:81
tuple tracks
Definition: testEve_cfg.py:39
HIPixelMedianVtxProducer(const edm::ParameterSet &ps)