CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
HIPixelMedianVtxProducer Class Reference

#include <HIPixelMedianVtxProducer.h>

Inheritance diagram for HIPixelMedianVtxProducer:
edm::one::EDProducer<> edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 HIPixelMedianVtxProducer (const edm::ParameterSet &ps)
 
void produce (edm::Event &ev, const edm::EventSetup &es) override
 
 ~HIPixelMedianVtxProducer () override
 
- Public Member Functions inherited from edm::one::EDProducer<>
 EDProducer ()=default
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
virtual ~ProducerBase () noexcept(false)
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void beginJob () override
 

Private Attributes

int theFitBinning
 
double theFitMaxZ
 
int theFitThreshold
 
int thePeakFindBinning
 
double thePeakFindMaxZ
 
unsigned int thePeakFindThresh
 
double thePtMin
 
edm::EDGetTokenT< reco::TrackCollectiontheTrackCollection
 

Additional Inherited Members

- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 12 of file HIPixelMedianVtxProducer.h.

Constructor & Destructor Documentation

HIPixelMedianVtxProducer::HIPixelMedianVtxProducer ( const edm::ParameterSet ps)
explicit

Definition at line 22 of file HIPixelMedianVtxProducer.cc.

References produce().

22  :
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 {
32  produces<reco::VertexCollection>();
33 
34  //In order to make fitting ROOT histograms thread safe
35  // one must call this undocumented function
36  TMinuitMinimizer::UseStaticMinuit(false);
37 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::TrackCollection > theTrackCollection
HIPixelMedianVtxProducer::~HIPixelMedianVtxProducer ( )
inlineoverride

Definition at line 16 of file HIPixelMedianVtxProducer.h.

References ev.

16 {};

Member Function Documentation

void HIPixelMedianVtxProducer::beginJob ( void  )
inlineoverrideprivatevirtual

Reimplemented from edm::one::EDProducerBase.

Definition at line 20 of file HIPixelMedianVtxProducer.h.

20 {};
void HIPixelMedianVtxProducer::produce ( edm::Event ev,
const edm::EventSetup es 
)
override

Definition at line 41 of file HIPixelMedianVtxProducer.cc.

References connectstrParser::f1, edm::Event::getByToken(), trackerHits::histo, mps_fire::i, LogTrace, eostools::move(), pileupDistInMC::num, edm::Handle< T >::product(), EnergyCorrector::pt, edm::Event::put(), theFitBinning, theFitMaxZ, theFitThreshold, thePeakFindBinning, thePeakFindMaxZ, thePeakFindThresh, thePtMin, theTrackCollection, HiIsolationCommonParameters_cff::track, findElectronsInSiStrips_cfi::trackCollection, l1t::tracks, tracks_, and electrons_cff::vertices.

Referenced by HIPixelMedianVtxProducer().

42 {
43  // Get pixel tracks
45  ev.getByToken(theTrackCollection, trackCollection);
46  const reco::TrackCollection tracks_ = *(trackCollection.product());
47 
48  // Select tracks above minimum pt
49  std::vector<const reco::Track *> tracks;
50  for (unsigned int i=0; i<tracks_.size(); i++) {
51  if (tracks_[i].pt() < thePtMin && std::fabs(tracks_[i].vz()) < 100000.) continue;
52  reco::TrackRef recTrack(trackCollection, i);
53  tracks.push_back( &(*recTrack));
54  }
55 
56  LogTrace("MinBiasTracking")
57  << " [VertexProducer] selected tracks: "
58  << tracks.size() << " (out of " << tracks_.size()
59  << ") above pt = " << thePtMin;
60 
61  // Output vertex collection
62  auto vertices = std::make_unique<reco::VertexCollection>();
63 
64  // No tracks -> return empty collection
65  if(tracks.empty()) {
66  ev.put(std::move(vertices));
67  return;
68  }
69 
70  // Sort tracks according to vertex z position
71  std::sort(tracks.begin(), tracks.end(), ComparePairs());
72 
73  // Calculate median vz
74  float med;
75  if(tracks.size() % 2 == 0)
76  med = (tracks[tracks.size()/2-1]->vz() + tracks[tracks.size()/2]->vz())/2;
77  else
78  med = tracks[tracks.size()/2 ]->vz();
79 
80  LogTrace("MinBiasTracking")
81  << " [vertex position] median = " << med << " cm";
82 
83  // In high multiplicity events, fit around most probable position
84  if(tracks.size() > thePeakFindThresh) {
85 
86  // Find maximum bin
87  TH1F hmax("hmax","hmax",thePeakFindBinning*2.0*thePeakFindMaxZ,-1.0*thePeakFindMaxZ,thePeakFindMaxZ);
88 
89  for(std::vector<const reco::Track *>::const_iterator
90  track = tracks.begin(); track!= tracks.end(); track++)
91  if(fabs((*track)->vz()) < thePeakFindMaxZ)
92  hmax.Fill((*track)->vz());
93 
94  int maxBin = hmax.GetMaximumBin();
95 
96  LogTrace("MinBiasTracking")
97  << " [vertex position] most prob = "
98  << hmax.GetBinCenter(maxBin) << " cm";
99 
100  // Find 3-bin weighted average
101  float num=0.0, denom=0.0;
102  for(int i=-1; i<=1; i++) {
103  num += hmax.GetBinContent(maxBin+i)*hmax.GetBinCenter(maxBin+i);
104  denom += hmax.GetBinContent(maxBin+i);
105  }
106 
107  if(denom==0)
108  {
110  err(2,2) = 0.1 * 0.1;
111  reco::Vertex ver(reco::Vertex::Point(0,0,99999.),
112  err, 0, 0, 1);
113  vertices->push_back(ver);
114  ev.put(std::move(vertices));
115  return;
116  }
117 
118  float nBinAvg = num/denom;
119 
120  // Center fit at 3-bin weighted average around max bin
121  med = nBinAvg;
122 
123  LogTrace("MinBiasTracking")
124  << " [vertex position] 3-bin weighted average = "
125  << nBinAvg << " cm";
126  }
127 
128  // Bin vz-values around most probable value (or median) for fitting
129  TH1F histo("histo","histo", theFitBinning*2.0*theFitMaxZ,-1.0*theFitMaxZ,theFitMaxZ);
130  histo.Sumw2();
131 
132  for(std::vector<const reco::Track *>::const_iterator
133  track = tracks.begin(); track!= tracks.end(); track++)
134  if(fabs((*track)->vz() - med) < theFitMaxZ)
135  histo.Fill((*track)->vz() - med);
136 
137  LogTrace("MinBiasTracking")
138  << " [vertex position] most prob for fit = "
139  << med + histo.GetBinCenter(histo.GetMaximumBin()) << " cm";
140 
141  // If there are very few entries, don't do the fit
142  if(histo.GetEntries() <= theFitThreshold) {
143 
144  LogTrace("MinBiasTracking")
145  << " [vertex position] Fewer than" << theFitThreshold
146  << " entries in fit histogram. Returning median.";
147 
149  err(2,2) = 0.1 * 0.1;
150  reco::Vertex ver(reco::Vertex::Point(0,0,med),
151  err, 0, 1, 1);
152  vertices->push_back(ver);
153  ev.put(std::move(vertices));
154  return;
155  }
156 
157  // Otherwise, there are enough entries to refine the estimate with a fit
158  TF1 f1("f1","[0]*exp(-0.5 * ((x-[1])/[2])^2) + [3]");
159  f1.SetParameters(10.,0.,0.02,0.002*tracks.size());
160  f1.SetParLimits(1,-0.1,0.1);
161  f1.SetParLimits(2,0.001,0.05);
162  f1.SetParLimits(3,0.0,0.005*tracks.size());
163 
164  histo.Fit(&f1,"QN");
165 
166  LogTrace("MinBiasTracking")
167  << " [vertex position] fitted = "
168  << med + f1.GetParameter(1) << " +- " << f1.GetParError(1) << " cm";
169 
171  err(2,2) = f1.GetParError(1) * f1.GetParError(1);
172  reco::Vertex ver(reco::Vertex::Point(0,0,med + f1.GetParameter(1)),
173  err, 0, 1, 1);
174  vertices->push_back(ver);
175 
176  ev.put(std::move(vertices));
177  return;
178 
179 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
const std::vector< reco::PFCandidatePtr > & tracks_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
#define LogTrace(id)
T const * product() const
Definition: Handle.h:81
def move(src, dest)
Definition: eostools.py:510
edm::EDGetTokenT< reco::TrackCollection > theTrackCollection

Member Data Documentation

int HIPixelMedianVtxProducer::theFitBinning
private

Definition at line 29 of file HIPixelMedianVtxProducer.h.

Referenced by produce().

double HIPixelMedianVtxProducer::theFitMaxZ
private

Definition at line 28 of file HIPixelMedianVtxProducer.h.

Referenced by produce().

int HIPixelMedianVtxProducer::theFitThreshold
private

Definition at line 27 of file HIPixelMedianVtxProducer.h.

Referenced by produce().

int HIPixelMedianVtxProducer::thePeakFindBinning
private

Definition at line 26 of file HIPixelMedianVtxProducer.h.

Referenced by produce().

double HIPixelMedianVtxProducer::thePeakFindMaxZ
private

Definition at line 25 of file HIPixelMedianVtxProducer.h.

Referenced by produce().

unsigned int HIPixelMedianVtxProducer::thePeakFindThresh
private

Definition at line 24 of file HIPixelMedianVtxProducer.h.

Referenced by produce().

double HIPixelMedianVtxProducer::thePtMin
private

Definition at line 23 of file HIPixelMedianVtxProducer.h.

Referenced by produce().

edm::EDGetTokenT<reco::TrackCollection> HIPixelMedianVtxProducer::theTrackCollection
private

Definition at line 20 of file HIPixelMedianVtxProducer.h.

Referenced by produce().