CMS 3D CMS Logo

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 #include "TMinuitMinimizer.h"
20 
21 /*****************************************************************************/
23  : theTrackCollection(consumes<reco::TrackCollection>(ps.getParameter<edm::InputTag>("TrackCollection"))),
24  thePtMin(ps.getParameter<double>("PtMin")),
25  thePeakFindThresh(ps.getParameter<unsigned int>("PeakFindThreshold")),
26  thePeakFindMaxZ(ps.getParameter<double>("PeakFindMaxZ")),
27  thePeakFindBinning(ps.getParameter<int>("PeakFindBinsPerCm")),
28  theFitThreshold(ps.getParameter<int>("FitThreshold")),
29  theFitMaxZ(ps.getParameter<double>("FitMaxZ")),
30  theFitBinning(ps.getParameter<int>("FitBinsPerCm")) {
31  produces<reco::VertexCollection>();
32 
33  //In order to make fitting ROOT histograms thread safe
34  // one must call this undocumented function
35  TMinuitMinimizer::UseStaticMinuit(false);
36 }
37 
38 /*****************************************************************************/
40  // Get pixel tracks
43  const reco::TrackCollection tracks_ = *(trackCollection.product());
44 
45  // Select tracks above minimum pt
46  std::vector<const reco::Track*> tracks;
47  for (unsigned int i = 0; i < tracks_.size(); i++) {
48  if (tracks_[i].pt() < thePtMin && std::fabs(tracks_[i].vz()) < 100000.)
49  continue;
50  reco::TrackRef recTrack(trackCollection, i);
51  tracks.push_back(&(*recTrack));
52  }
53 
54  LogTrace("MinBiasTracking") << " [VertexProducer] selected tracks: " << tracks.size() << " (out of " << tracks_.size()
55  << ") above pt = " << thePtMin;
56 
57  // Output vertex collection
58  auto vertices = std::make_unique<reco::VertexCollection>();
59 
60  // No tracks -> return empty collection
61  if (tracks.empty()) {
62  ev.put(std::move(vertices));
63  return;
64  }
65 
66  // Sort tracks according to vertex z position
67  std::sort(tracks.begin(), tracks.end(), ComparePairs());
68 
69  // Calculate median vz
70  float med;
71  if (tracks.size() % 2 == 0)
72  med = (tracks[tracks.size() / 2 - 1]->vz() + tracks[tracks.size() / 2]->vz()) / 2;
73  else
74  med = tracks[tracks.size() / 2]->vz();
75 
76  LogTrace("MinBiasTracking") << " [vertex position] median = " << med << " cm";
77 
78  // In high multiplicity events, fit around most probable position
79  if (tracks.size() > thePeakFindThresh) {
80  // Find maximum bin
81  TH1F hmax("hmax", "hmax", thePeakFindBinning * 2.0 * thePeakFindMaxZ, -1.0 * thePeakFindMaxZ, thePeakFindMaxZ);
82 
83  for (std::vector<const reco::Track*>::const_iterator track = tracks.begin(); track != tracks.end(); track++)
84  if (fabs((*track)->vz()) < thePeakFindMaxZ)
85  hmax.Fill((*track)->vz());
86 
87  int maxBin = hmax.GetMaximumBin();
88 
89  LogTrace("MinBiasTracking") << " [vertex position] most prob = " << hmax.GetBinCenter(maxBin) << " cm";
90 
91  // Find 3-bin weighted average
92  float num = 0.0, denom = 0.0;
93  for (int i = -1; i <= 1; i++) {
94  num += hmax.GetBinContent(maxBin + i) * hmax.GetBinCenter(maxBin + i);
95  denom += hmax.GetBinContent(maxBin + i);
96  }
97 
98  if (denom == 0) {
100  err(2, 2) = 0.1 * 0.1;
101  reco::Vertex ver(reco::Vertex::Point(0, 0, 99999.), err, 0, 0, 1);
102  vertices->push_back(ver);
103  ev.put(std::move(vertices));
104  return;
105  }
106 
107  float nBinAvg = num / denom;
108 
109  // Center fit at 3-bin weighted average around max bin
110  med = nBinAvg;
111 
112  LogTrace("MinBiasTracking") << " [vertex position] 3-bin weighted average = " << nBinAvg << " cm";
113  }
114 
115  // Bin vz-values around most probable value (or median) for fitting
116  TH1F histo("histo", "histo", theFitBinning * 2.0 * theFitMaxZ, -1.0 * theFitMaxZ, theFitMaxZ);
117  histo.Sumw2();
118 
119  for (std::vector<const reco::Track*>::const_iterator track = tracks.begin(); track != tracks.end(); track++)
120  if (fabs((*track)->vz() - med) < theFitMaxZ)
121  histo.Fill((*track)->vz() - med);
122 
123  LogTrace("MinBiasTracking") << " [vertex position] most prob for fit = "
124  << med + histo.GetBinCenter(histo.GetMaximumBin()) << " cm";
125 
126  // If there are very few entries, don't do the fit
127  if (histo.GetEntries() <= theFitThreshold) {
128  LogTrace("MinBiasTracking") << " [vertex position] Fewer than" << theFitThreshold
129  << " entries in fit histogram. Returning median.";
130 
132  err(2, 2) = 0.1 * 0.1;
133  reco::Vertex ver(reco::Vertex::Point(0, 0, med), err, 0, 1, 1);
134  vertices->push_back(ver);
135  ev.put(std::move(vertices));
136  return;
137  }
138 
139  // Otherwise, there are enough entries to refine the estimate with a fit
140  TF1 f1("f1", "[0]*exp(-0.5 * ((x-[1])/[2])^2) + [3]");
141  f1.SetParameters(10., 0., 0.02, 0.002 * tracks.size());
142  f1.SetParLimits(1, -0.1, 0.1);
143  f1.SetParLimits(2, 0.001, 0.05);
144  f1.SetParLimits(3, 0.0, 0.005 * tracks.size());
145 
146  histo.Fit(&f1, "QN SERIAL");
147 
148  LogTrace("MinBiasTracking") << " [vertex position] fitted = " << med + f1.GetParameter(1) << " +- "
149  << f1.GetParError(1) << " cm";
150 
152  err(2, 2) = f1.GetParError(1) * f1.GetParError(1);
153  reco::Vertex ver(reco::Vertex::Point(0, 0, med + f1.GetParameter(1)), err, 0, 1, 1);
154  vertices->push_back(ver);
155 
156  ev.put(std::move(vertices));
157  return;
158 }
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
mps_fire.i
i
Definition: mps_fire.py:355
HIPixelMedianVtxProducer::theFitThreshold
int theFitThreshold
Definition: HIPixelMedianVtxProducer.h:29
makePileupJSON.denom
denom
Definition: makePileupJSON.py:147
ESHandle.h
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
edm
HLT enums.
Definition: AlignableModifier.h:19
reco::Vertex::Error
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
HIPixelMedianVtxProducer::thePeakFindBinning
int thePeakFindBinning
Definition: HIPixelMedianVtxProducer.h:28
timingPdfMaker.histo
histo
Definition: timingPdfMaker.py:279
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
HIPixelMedianVtxProducer::thePeakFindMaxZ
double thePeakFindMaxZ
Definition: HIPixelMedianVtxProducer.h:27
tools.TF1
TF1
Definition: tools.py:23
edm::Handle< reco::TrackCollection >
edm::Ref< TrackCollection >
MakerMacros.h
HIPixelMedianVtxProducer::theTrackCollection
edm::EDGetTokenT< reco::TrackCollection > theTrackCollection
Definition: HIPixelMedianVtxProducer.h:22
ComparePairs
Definition: HIPixelMedianVtxProducer.h:34
Vertex.h
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
HIPixelMedianVtxProducer::theFitMaxZ
double theFitMaxZ
Definition: HIPixelMedianVtxProducer.h:30
duplicaterechits_cfi.trackCollection
trackCollection
Definition: duplicaterechits_cfi.py:4
Event.h
HIPixelMedianVtxProducer::HIPixelMedianVtxProducer
HIPixelMedianVtxProducer(const edm::ParameterSet &ps)
Definition: HIPixelMedianVtxProducer.cc:22
runTheMatrix.err
err
Definition: runTheMatrix.py:288
cms::cuda::allocator::maxBin
constexpr unsigned int maxBin
Definition: getCachingDeviceAllocator.h:20
createfilelist.int
int
Definition: createfilelist.py:10
HIPixelMedianVtxProducer::produce
void produce(edm::Event &ev, const edm::EventSetup &es) override
Definition: HIPixelMedianVtxProducer.cc:39
edm::EventSetup
Definition: EventSetup.h:57
HIPixelMedianVtxProducer::thePeakFindThresh
unsigned int thePeakFindThresh
Definition: HIPixelMedianVtxProducer.h:26
EgammaValidation_cff.num
num
Definition: EgammaValidation_cff.py:34
VertexFwd.h
HIPixelMedianVtxProducer::thePtMin
double thePtMin
Definition: HIPixelMedianVtxProducer.h:25
reco::Vertex::Point
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
eostools.move
def move(src, dest)
Definition: eostools.py:511
Frameworkfwd.h
HIPixelMedianVtxProducer::theFitBinning
int theFitBinning
Definition: HIPixelMedianVtxProducer.h:31
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
LogTrace
#define LogTrace(id)
Definition: MessageLogger.h:671
edm::Event
Definition: Event.h:73
DeadROC_duringRun.f1
f1
Definition: DeadROC_duringRun.py:219
HIPixelMedianVtxProducer.h
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
reco::Vertex
Definition: Vertex.h:35
pwdgSkimBPark_cfi.vertices
vertices
Definition: pwdgSkimBPark_cfi.py:7