CMS 3D CMS Logo

PixelVertexProducerMedian.cc
Go to the documentation of this file.
2 
7 
10 
13 
14 #include <fstream>
15 #include <iostream>
16 #include <vector>
17 #include <algorithm>
18 
19 #include "TROOT.h"
20 #include "TH1F.h"
21 #include "TF1.h"
22 
23 /*****************************************************************************/
24 struct ComparePairs {
25  bool operator()(const reco::Track* t1, const reco::Track* t2) { return (t1->vz() < t2->vz()); };
26 };
27 
28 /*****************************************************************************/
30  thePtMin = theConfig.getParameter<double>("PtMin");
31  produces<reco::VertexCollection>();
32 }
33 
34 /*****************************************************************************/
36 
37 /*****************************************************************************/
39  // Get pixel tracks
41  std::string trackCollectionName = theConfig.getParameter<std::string>("TrackCollection");
42  ev.getByLabel(trackCollectionName, trackCollection);
43  const reco::TrackCollection tracks_ = *(trackCollection.product());
44 
45  // Select tracks
46  std::vector<const reco::Track*> tracks;
47  for (unsigned int i = 0; i < tracks_.size(); i++) {
48  if (tracks_[i].pt() > thePtMin) {
49  reco::TrackRef recTrack(trackCollection, i);
50  tracks.push_back(&(*recTrack));
51  }
52  }
53 
54  LogTrace("MinBiasTracking") << " [VertexProducer] selected tracks: " << tracks.size() << " (out of " << tracks_.size()
55  << ")";
56 
57  auto vertices = std::make_unique<reco::VertexCollection>();
58 
59  if (!tracks.empty()) {
60  // Sort along vertex z position
61  std::sort(tracks.begin(), tracks.end(), ComparePairs());
62 
63  // Median
64  float med;
65  if (tracks.size() % 2 == 0)
66  med = (tracks[tracks.size() / 2 - 1]->vz() + tracks[tracks.size() / 2]->vz()) / 2;
67  else
68  med = tracks[tracks.size() / 2]->vz();
69 
70  LogTrace("MinBiasTracking") << " [vertex position] median = " << med << " cm";
71 
72  if (tracks.size() > 10) {
73  // Binning around med, halfWidth
74  int nBin = 100;
75  float halfWidth = 0.1; // cm
76 
77  // Most probable
78  TH1F histo("histo", "histo", nBin, -halfWidth, halfWidth);
79 
80  for (std::vector<const reco::Track*>::const_iterator track = tracks.begin(); track != tracks.end(); track++)
81  if (fabs((*track)->vz() - med) < halfWidth)
82  histo.Fill((*track)->vz() - med);
83 
84  LogTrace("MinBiasTracking") << " [vertex position] most prob = "
85  << med + histo.GetBinCenter(histo.GetMaximumBin()) << " cm";
86 
87  // Fit above max/2
88  histo.Sumw2();
89 
90  TF1 f1("f1", "[0]*exp(-0.5 * ((x-[1])/[2])^2) + [3]");
91  f1.SetParameters(10., 0., 0.01, 1.);
92 
93  histo.Fit(&f1, "QN");
94 
95  LogTrace("MinBiasTracking") << " [vertex position] fitted = " << med + f1.GetParameter(1) << " +- "
96  << f1.GetParError(1) << " cm";
97 
98  // Store
100  err(2, 2) = f1.GetParError(1) * f1.GetParError(1);
101  reco::Vertex ver(reco::Vertex::Point(0, 0, med + f1.GetParameter(1)), err, 0, 1, 1);
102  vertices->push_back(ver);
103  } else {
104  // Store
106  err(2, 2) = 0.1 * 0.1;
107  reco::Vertex ver(reco::Vertex::Point(0, 0, med), err, 0, 1, 1);
108  vertices->push_back(ver);
109  }
110  }
111  ev.put(std::move(vertices));
112 }
edm::StreamID
Definition: StreamID.h:30
RandomServiceHelper.t2
t2
Definition: RandomServiceHelper.py:257
mps_fire.i
i
Definition: mps_fire.py:428
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11724
ESHandle.h
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
PixelVertexProducerMedian::~PixelVertexProducerMedian
~PixelVertexProducerMedian() override
Definition: PixelVertexProducerMedian.cc:35
reco::Vertex::Error
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
timingPdfMaker.histo
histo
Definition: timingPdfMaker.py:278
PixelVertexProducerMedian.h
edm::Handle< reco::TrackCollection >
PFElectronDQMAnalyzer_cfi.nBin
nBin
Definition: PFElectronDQMAnalyzer_cfi.py:25
AlignmentTracksFromVertexSelector_cfi.vertices
vertices
Definition: AlignmentTracksFromVertexSelector_cfi.py:5
edm::Ref< TrackCollection >
PixelVertexProducerMedian::thePtMin
double thePtMin
Definition: PixelVertexProducerMedian.h:20
MakerMacros.h
RandomServiceHelper.t1
t1
Definition: RandomServiceHelper.py:256
Track.h
TrackFwd.h
PixelVertexProducerMedian::produce
void produce(edm::StreamID, edm::Event &ev, const edm::EventSetup &es) const override
Definition: PixelVertexProducerMedian.cc:38
reco::Track
Definition: Track.h:27
ComparePairs
Definition: HIPixelMedianVtxProducer.h:34
edm::ParameterSet
Definition: ParameterSet.h:47
duplicaterechits_cfi.trackCollection
trackCollection
Definition: duplicaterechits_cfi.py:4
Event.h
tracks
const uint32_t *__restrict__ const HitContainer *__restrict__ TkSoA *__restrict__ tracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:176
jetUpdater_cfi.sort
sort
Definition: jetUpdater_cfi.py:29
edm::EventSetup
Definition: EventSetup.h:58
submitPVResolutionJobs.err
err
Definition: submitPVResolutionJobs.py:85
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
PixelVertexProducerMedian::PixelVertexProducerMedian
PixelVertexProducerMedian(const edm::ParameterSet &ps)
Definition: PixelVertexProducerMedian.cc:29
VertexFwd.h
reco::Vertex::Point
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
eostools.move
def move(src, dest)
Definition: eostools.py:511
Vertex.h
Frameworkfwd.h
ev
bool ev
Definition: Hydjet2Hadronizer.cc:97
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
LogTrace
#define LogTrace(id)
Definition: MessageLogger.h:234
ComparePairs::operator()
bool operator()(const reco::Track *t1, const reco::Track *t2)
Definition: PixelVertexProducerMedian.cc:25
edm::Event
Definition: Event.h:73
DeadROC_duringRun.f1
f1
Definition: DeadROC_duringRun.py:219
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
reco::Vertex
Definition: Vertex.h:35
PixelVertexProducerMedian::theConfig
edm::ParameterSet theConfig
Definition: PixelVertexProducerMedian.h:19