CMS 3D CMS Logo

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

#include <RecoBTag/FastPrimaryVertexWithWeightsProducer/src/FastPrimaryVertexWithWeightsProducer.cc>

Inheritance diagram for FastPrimaryVertexWithWeightsProducer:
edm::stream::EDProducer<>

Public Member Functions

 FastPrimaryVertexWithWeightsProducer (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

edm::EDGetTokenT< reco::BeamSpotbeamSpotToken
 
edm::EDGetTokenT< SiPixelClusterCollectionNewclustersToken
 
edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken
 
const bool m_barrel
 
const edm::InputTag m_clusters
 
const double m_EC_weight
 
const bool m_endCap
 
const double m_maxDeltaPhi
 
const double m_maxDeltaPhi_EC
 
const double m_maxJetEta
 
const double m_maxJetEta_EC
 
const double m_maxSizeX
 
const double m_maxSizeY_q
 
const double m_maxZ
 
const double m_minJetEta_EC
 
const double m_minJetPt
 
const double m_minSizeY_q
 
const int m_njets
 
const double m_peakSizeY_q
 
const double m_PixelCellHeightOverWidth
 
std::string m_pixelCPE
 
const bool m_ptWeighting
 
const double m_ptWeighting_offset
 
const double m_ptWeighting_slope
 
const double m_weight_charge_down
 
const double m_weight_charge_peak
 
const double m_weight_charge_up
 
const double m_weight_dPhi
 
const double m_weight_dPhi_EC
 
const double m_weight_rho_up
 
const double m_weight_SizeX1
 
const double m_weightCut_step2
 
const double m_weightCut_step3
 
const double m_zClusterSearchArea_step2
 
const double m_zClusterSearchArea_step3
 
const double m_zClusterWidth_step1
 
const double m_zClusterWidth_step2
 
const double m_zClusterWidth_step3
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Description: The FastPrimaryVertex is an algorithm used to find the primary vertex of an event @HLT. It takes the pixel clusters compabible with a jet and project it to the beamSpot with the eta-angle of the jet. The peak on the z-projected clusters distribution is our FastPrimaryVertex.

Implementation: [Notes on implementation]

Definition at line 56 of file FastPrimaryVertexWithWeightsProducer.cc.

Constructor & Destructor Documentation

◆ FastPrimaryVertexWithWeightsProducer()

FastPrimaryVertexWithWeightsProducer::FastPrimaryVertexWithWeightsProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 127 of file FastPrimaryVertexWithWeightsProducer.cc.

128  : m_maxZ(iConfig.getParameter<double>("maxZ")),
129  m_pixelCPE(iConfig.getParameter<std::string>("pixelCPE")),
130 
131  m_njets(iConfig.getParameter<int>("njets")),
132  m_maxJetEta(iConfig.getParameter<double>("maxJetEta")),
133  m_minJetPt(iConfig.getParameter<double>("minJetPt")),
134 
135  m_barrel(iConfig.getParameter<bool>("barrel")),
136  m_maxSizeX(iConfig.getParameter<double>("maxSizeX")),
137  m_maxDeltaPhi(iConfig.getParameter<double>("maxDeltaPhi")),
138  m_weight_charge_down(iConfig.getParameter<double>("weight_charge_down")),
139  m_weight_charge_up(iConfig.getParameter<double>("weight_charge_up")),
140  m_PixelCellHeightOverWidth(iConfig.getParameter<double>("PixelCellHeightOverWidth")),
141  m_minSizeY_q(iConfig.getParameter<double>("minSizeY_q")),
142  m_maxSizeY_q(iConfig.getParameter<double>("maxSizeY_q")),
143 
144  m_weight_dPhi(iConfig.getParameter<double>("weight_dPhi")),
145  m_weight_SizeX1(iConfig.getParameter<double>("weight_SizeX1")),
146  m_weight_rho_up(iConfig.getParameter<double>("weight_rho_up")),
147  m_weight_charge_peak(iConfig.getParameter<double>("weight_charge_peak")),
148  m_peakSizeY_q(iConfig.getParameter<double>("peakSizeY_q")),
149 
150  m_endCap(iConfig.getParameter<bool>("endCap")),
151  m_minJetEta_EC(iConfig.getParameter<double>("minJetEta_EC")),
152  m_maxJetEta_EC(iConfig.getParameter<double>("maxJetEta_EC")),
153  m_maxDeltaPhi_EC(iConfig.getParameter<double>("maxDeltaPhi_EC")),
154  m_EC_weight(iConfig.getParameter<double>("EC_weight")),
155  m_weight_dPhi_EC(iConfig.getParameter<double>("weight_dPhi_EC")),
156 
157  m_zClusterWidth_step1(iConfig.getParameter<double>("zClusterWidth_step1")),
158 
159  m_zClusterWidth_step2(iConfig.getParameter<double>("zClusterWidth_step2")),
160  m_zClusterSearchArea_step2(iConfig.getParameter<double>("zClusterSearchArea_step2")),
161  m_weightCut_step2(iConfig.getParameter<double>("weightCut_step2")),
162 
163  m_zClusterWidth_step3(iConfig.getParameter<double>("zClusterWidth_step3")),
164  m_zClusterSearchArea_step3(iConfig.getParameter<double>("zClusterSearchArea_step3")),
165  m_weightCut_step3(iConfig.getParameter<double>("weightCut_step3")),
166 
167  m_ptWeighting(iConfig.getParameter<bool>("ptWeighting")),
168  m_ptWeighting_slope(iConfig.getParameter<double>("ptWeighting_slope")),
169  m_ptWeighting_offset(iConfig.getParameter<double>("ptWeighting_offset")) {
170  clustersToken = consumes<SiPixelClusterCollectionNew>(iConfig.getParameter<edm::InputTag>("clusters"));
171  beamSpotToken = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"));
172  jetsToken = consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("jets"));
173 
174  produces<reco::VertexCollection>();
175  produces<float>();
176 }

References beamSpotToken, clustersToken, edm::ParameterSet::getParameter(), and jetsToken.

Member Function Documentation

◆ fillDescriptions()

void FastPrimaryVertexWithWeightsProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 178 of file FastPrimaryVertexWithWeightsProducer.cc.

178  {
180  desc.add<edm::InputTag>("clusters", edm::InputTag("hltSiPixelClusters"));
181  desc.add<edm::InputTag>("beamSpot", edm::InputTag("hltOnlineBeamSpot"));
182  desc.add<edm::InputTag>("jets", edm::InputTag("hltCaloJetL1FastJetCorrected"));
183  desc.add<std::string>("pixelCPE", "hltESPPixelCPEGeneric");
184  desc.add<double>("maxZ", 19.0);
185  desc.add<int>("njets", 999);
186  desc.add<double>("maxJetEta", 2.6);
187  desc.add<double>("minJetPt", 40.);
188  desc.add<bool>("barrel", true);
189  desc.add<double>("maxSizeX", 2.1);
190  desc.add<double>("maxDeltaPhi", 0.21);
191  desc.add<double>("PixelCellHeightOverWidth", 1.8);
192  desc.add<double>("weight_charge_down", 11. * 1000.);
193  desc.add<double>("weight_charge_up", 190. * 1000.);
194  desc.add<double>("maxSizeY_q", 2.0);
195  desc.add<double>("minSizeY_q", -0.6);
196  desc.add<double>("weight_dPhi", 0.13888888);
197  desc.add<double>("weight_SizeX1", 0.88);
198  desc.add<double>("weight_rho_up", 22.);
199  desc.add<double>("weight_charge_peak", 22. * 1000.);
200  desc.add<double>("peakSizeY_q", 1.0);
201  desc.add<bool>("endCap", true);
202  desc.add<double>("minJetEta_EC", 1.3);
203  desc.add<double>("maxJetEta_EC", 2.6);
204  desc.add<double>("maxDeltaPhi_EC", 0.14);
205  desc.add<double>("EC_weight", 0.008);
206  desc.add<double>("weight_dPhi_EC", 0.064516129);
207  desc.add<double>("zClusterWidth_step1", 2.0);
208  desc.add<double>("zClusterWidth_step2", 0.65);
209  desc.add<double>("zClusterSearchArea_step2", 3.0);
210  desc.add<double>("weightCut_step2", 0.05);
211  desc.add<double>("zClusterWidth_step3", 0.3);
212  desc.add<double>("zClusterSearchArea_step3", 0.55);
213  desc.add<double>("weightCut_step3", 0.1);
214  desc.add<bool>("ptWeighting", false); // <---- newMethod
215  desc.add<double>("ptWeighting_slope", 1 / 20.);
216  desc.add<double>("ptWeighting_offset", -1);
217  descriptions.add("fastPrimaryVertexWithWeightsProducer", desc);
218 }

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), HLT_2018_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

◆ produce()

void FastPrimaryVertexWithWeightsProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

if it is a cluster to project

Definition at line 220 of file FastPrimaryVertexWithWeightsProducer.cc.

220  {
221  using namespace edm;
222  using namespace reco;
223  using namespace std;
224 
225  const float barrel_lenght = 30; //half lenght of the pixel barrel 30 cm
226 
227  //get pixel cluster
229  iEvent.getByToken(clustersToken, cH);
231 
232  //get jets
234  iEvent.getByToken(jetsToken, jH);
235  const edm::View<reco::Jet>& jets = *jH.product();
236 
237  vector<const reco::Jet*> selectedJets;
238  int countjet = 0;
239  for (edm::View<reco::Jet>::const_iterator it = jets.begin(); it != jets.end() && countjet < m_njets; it++) {
240  if ( //select jets used in barrel or endcap pixel cluster projections
241  ((it->pt() >= m_minJetPt) && std::abs(it->eta()) <= m_maxJetEta) || //barrel
242  ((it->pt() >= m_minJetPt) && std::abs(it->eta()) <= m_maxJetEta_EC &&
243  std::abs(it->eta()) >= m_minJetEta_EC) //endcap
244  ) {
245  selectedJets.push_back(&(*it));
246  countjet++;
247  }
248  }
249 
250  //get PixelClusterParameterEstimator
253  iSetup.get<TkPixelCPERecord>().get(m_pixelCPE, pe);
254  pp = pe.product();
255 
256  //get beamSpot
258  iEvent.getByToken(beamSpotToken, beamSpot);
259 
260  //get TrackerGeometry
263  const TrackerGeometry* trackerGeometry = tracker.product();
264 
265  // PART I: get z-projections with z-weights
266  std::vector<float> zProjections;
267  std::vector<float> zWeights;
268  int jet_count = 0;
269  for (vector<const reco::Jet*>::iterator jit = selectedJets.begin(); jit != selectedJets.end();
270  jit++) { //loop on selected jets
271  float px = (*jit)->px();
272  float py = (*jit)->py();
273  float pz = (*jit)->pz();
274  float pt = (*jit)->pt();
275  float eta = (*jit)->eta();
276  float jetZOverRho = (*jit)->momentum().Z() / (*jit)->momentum().Rho();
277  float pt_weight = pt * m_ptWeighting_slope + m_ptWeighting_offset;
279  it++) //Loop on pixel modules with clusters
280  { //loop on pixel modules
281  DetId id = it->detId();
282  const edmNew::DetSet<SiPixelCluster>& detset = (*it);
283  Point3DBase<float, GlobalTag> modulepos = trackerGeometry->idToDet(id)->position();
284  float zmodule = modulepos.z() -
285  ((modulepos.x() - beamSpot->x0()) * px + (modulepos.y() - beamSpot->y0()) * py) / pt * pz / pt;
286  if ((std::abs(deltaPhi((*jit)->momentum().Phi(), modulepos.phi())) < m_maxDeltaPhi * 2) &&
287  (std::abs(zmodule) < (m_maxZ + barrel_lenght))) { //if it is a compatible module
288  for (size_t j = 0; j < detset.size(); j++) { //loop on pixel clusters on this module
289  const SiPixelCluster& aCluster = detset[j];
290  if ( //it is a cluster to project
291  ( // barrel
292  m_barrel && std::abs(modulepos.z()) < barrel_lenght && pt >= m_minJetPt && jet_count < m_njets &&
293  std::abs(eta) <= m_maxJetEta && aCluster.sizeX() <= m_maxSizeX &&
294  aCluster.sizeY() >= m_PixelCellHeightOverWidth * std::abs(jetZOverRho) + m_minSizeY_q &&
295  aCluster.sizeY() <= m_PixelCellHeightOverWidth * std::abs(jetZOverRho) + m_maxSizeY_q) ||
296  ( // EC
297  m_endCap && std::abs(modulepos.z()) > barrel_lenght && pt > m_minJetPt && jet_count < m_njets &&
299  aCluster.sizeX() <= m_maxSizeX)) {
300  Point3DBase<float, GlobalTag> v = trackerGeometry->idToDet(id)->surface().toGlobal(
301  pp->localParametersV(aCluster, (*trackerGeometry->idToDetUnit(id)))[0].first);
302  GlobalPoint v_bs(v.x() - beamSpot->x0(), v.y() - beamSpot->y0(), v.z());
303  if ( //it pass DeltaPhi(Jet,Cluster) requirements
304  (m_barrel && std::abs(modulepos.z()) < barrel_lenght &&
305  std::abs(deltaPhi((*jit)->momentum().Phi(), v_bs.phi())) <= m_maxDeltaPhi) || //barrel
306  (m_endCap && std::abs(modulepos.z()) > barrel_lenght &&
307  std::abs(deltaPhi((*jit)->momentum().Phi(), v_bs.phi())) <= m_maxDeltaPhi_EC) //EC
308  ) {
309  //calculate z-projection
310  float z = v.z() - ((v.x() - beamSpot->x0()) * px + (v.y() - beamSpot->y0()) * py) / pt * pz / pt;
311  if (std::abs(z) < m_maxZ) {
312  zProjections.push_back(z); //add z-projection in zProjections
313  float weight = 0;
314  //calculate zWeight
315  if (std::abs(modulepos.z()) < barrel_lenght) { //barrel
316  //calculate weight_sizeY
317  float sizeY = aCluster.sizeY();
318  float sizeY_up = m_PixelCellHeightOverWidth * std::abs(jetZOverRho) + m_maxSizeY_q;
319  float sizeY_peak = m_PixelCellHeightOverWidth * std::abs(jetZOverRho) + m_peakSizeY_q;
320  float sizeY_down = m_PixelCellHeightOverWidth * std::abs(jetZOverRho) + m_minSizeY_q;
321  float weight_sizeY_up = (sizeY_up - sizeY) / (sizeY_up - sizeY_peak);
322  float weight_sizeY_down = (sizeY - sizeY_down) / (sizeY_peak - sizeY_down);
323  weight_sizeY_down = weight_sizeY_down * (weight_sizeY_down > 0) * (weight_sizeY_down < 1);
324  weight_sizeY_up = weight_sizeY_up * (weight_sizeY_up > 0) * (weight_sizeY_up < 1);
325  float weight_sizeY = weight_sizeY_up + weight_sizeY_down;
326 
327  //calculate weight_rho
328  float rho = sqrt(v_bs.x() * v_bs.x() + v_bs.y() * v_bs.y());
329  float weight_rho = ((m_weight_rho_up - rho) / m_weight_rho_up);
330 
331  //calculate weight_dPhi
332  float weight_dPhi = exp(-std::abs(deltaPhi((*jit)->momentum().Phi(), v_bs.phi())) / m_weight_dPhi);
333 
334  //calculate weight_sizeX1
335  float weight_sizeX1 = (aCluster.sizeX() == 2) + (aCluster.sizeX() == 1) * m_weight_SizeX1;
336 
337  //calculate weight_charge
338  float charge = aCluster.charge();
339  float weightCluster_up = (m_weight_charge_up - charge) / (m_weight_charge_up - m_weight_charge_peak);
340  float weightCluster_down =
342  weightCluster_down = weightCluster_down * (weightCluster_down > 0) * (weightCluster_down < 1);
343  weightCluster_up = weightCluster_up * (weightCluster_up > 0) * (weightCluster_up < 1);
344  float weight_charge = weightCluster_up + weightCluster_down;
345 
346  //calculate the final weight
347  weight = weight_dPhi * weight_sizeY * weight_rho * weight_sizeX1 * weight_charge;
348  } else if (std::abs(modulepos.z()) > barrel_lenght) // EC
349  { // EC
350  //calculate weight_dPhi
351  float weight_dPhi = exp(-std::abs(deltaPhi((*jit)->momentum().Phi(), v_bs.phi())) / m_weight_dPhi_EC);
352  //calculate the final weight
354  }
355  if (m_ptWeighting)
356  weight = weight * pt_weight;
357  zWeights.push_back(weight); //add the weight to zWeights
358  }
359  } //if it pass DeltaPhi(Jet,Cluster) requirements
360  }
361  } //loop on pixel clusters on this module
362  } //if it is a compatible module
363  } //loop on pixel modules
364  jet_count++;
365  } //loop on selected jets
366 
367  //order zProjections and zWeights by z
368  std::multimap<float, float> zWithW;
369  size_t i = 0;
370  for (i = 0; i < zProjections.size(); i++)
371  zWithW.insert(std::pair<float, float>(zProjections[i], zWeights[i]));
372  i = 0;
373  for (std::multimap<float, float>::iterator it = zWithW.begin(); it != zWithW.end(); it++, i++) {
374  zProjections[i] = it->first;
375  zWeights[i] = it->second;
376  } //order zProjections and zWeights by z
377 
378  //calculate zWeightsSquared
379  std::vector<float> zWeightsSquared;
380  for (std::vector<float>::iterator it = zWeights.begin(); it != zWeights.end(); it++) {
381  zWeightsSquared.push_back((*it) * (*it));
382  }
383 
384  //do multi-step peak searching
385  float res_step1 = FindPeakFastPV(zProjections, zWeights, 0.0, m_zClusterWidth_step1, 999.0, -1.0);
386  float res_step2 = FindPeakFastPV(
387  zProjections, zWeights, res_step1, m_zClusterWidth_step2, m_zClusterSearchArea_step2, m_weightCut_step2);
388  float res_step3 = FindPeakFastPV(zProjections,
389  zWeightsSquared,
390  res_step2,
394 
395  float centerWMax = res_step3;
396  //End of PART II
397 
398  //Make the output
399  float res = 0;
400  if (zProjections.size() > 2) {
401  res = centerWMax;
403  e(0, 0) = 0.0015 * 0.0015;
404  e(1, 1) = 0.0015 * 0.0015;
405  e(2, 2) = 1.5 * 1.5;
407  Vertex thePV(p, e, 1, 1, 0);
408  auto pOut = std::make_unique<reco::VertexCollection>();
409  pOut->push_back(thePV);
410  iEvent.put(std::move(pOut));
411  } else {
413  e(0, 0) = 0.0015 * 0.0015;
414  e(1, 1) = 0.0015 * 0.0015;
415  e(2, 2) = 1.5 * 1.5;
417  Vertex thePV(p, e, 0, 0, 0);
418  auto pOut = std::make_unique<reco::VertexCollection>();
419  pOut->push_back(thePV);
420  iEvent.put(std::move(pOut));
421  }
422 
423  //Finally, calculate the zClusterQuality as Sum(weights near the fastPV)/sqrt(Sum(total weights)) [a kind of zCluster significance]
424 
425  const float half_width_peak = 1;
426  float nWeightedTot = 0;
427  float nWeightedTotPeak = 0;
428  for (std::vector<float>::iterator it = zProjections.begin(); it != zProjections.end(); it++) {
429  nWeightedTot += zWeights[it - zProjections.begin()];
430  if ((res - half_width_peak) <= (*it) && (*it) <= (res + half_width_peak)) {
431  nWeightedTotPeak += zWeights[it - zProjections.begin()];
432  }
433  }
434 
435  auto zClusterQuality = std::make_unique<float>();
436  *zClusterQuality = -1;
437  if (nWeightedTot != 0) {
438  // where 30 is the beam spot lenght
439  *zClusterQuality = nWeightedTotPeak / sqrt(nWeightedTot / (2 * half_width_peak));
440  iEvent.put(std::move(zClusterQuality));
441  } else
442  iEvent.put(std::move(zClusterQuality));
443 }

References funct::abs(), pwdgSkimBPark_cfi::beamSpot, beamSpotToken, ALCARECOTkAlJpsiMuMu_cff::charge, SiPixelCluster::charge(), clustersToken, SiPixelRawToDigiRegional_cfi::deltaPhi, MillePedeFileConverter_cfg::e, PVValHelper::eta, JetChargeProducer_cfi::exp, FindPeakFastPV(), edm::EventSetup::get(), get, mps_fire::i, TrackerGeometry::idToDet(), TrackerGeometry::idToDetUnit(), iEvent, dqmiolumiharvest::j, singleTopDQM_cfi::jets, jetsToken, m_barrel, m_EC_weight, m_endCap, m_maxDeltaPhi, m_maxDeltaPhi_EC, m_maxJetEta, m_maxJetEta_EC, m_maxSizeX, m_maxSizeY_q, m_maxZ, m_minJetEta_EC, m_minJetPt, m_minSizeY_q, m_njets, m_peakSizeY_q, m_PixelCellHeightOverWidth, m_pixelCPE, m_ptWeighting, m_ptWeighting_offset, m_ptWeighting_slope, m_weight_charge_down, m_weight_charge_peak, m_weight_charge_up, m_weight_dPhi, m_weight_dPhi_EC, m_weight_rho_up, m_weight_SizeX1, m_weightCut_step2, m_weightCut_step3, m_zClusterSearchArea_step2, m_zClusterSearchArea_step3, m_zClusterWidth_step1, m_zClusterWidth_step2, m_zClusterWidth_step3, eostools::move(), AlCaHLTBitMon_ParallelJobs::p, PV3DBase< T, PVType, FrameType >::phi(), LumiMonitor_cff::pixelClusters, GeomDet::position(), createTree::pp, edm::Handle< T >::product(), edm::ESHandle< T >::product(), DiDispStaMuonMonitor_cfi::pt, multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, rho, objectSelection_cff::selectedJets, edmNew::DetSet< T >::size(), SiPixelCluster::sizeX(), SiPixelCluster::sizeY(), RecoTauValidation_cfi::sizeY, mathSSE::sqrt(), GeomDet::surface(), Surface::toGlobal(), PbPb_ZMuSkimMuonDPG_cff::tracker, findQualityFiles::v, HLT_2018_cff::weight_dPhi, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

Member Data Documentation

◆ beamSpotToken

edm::EDGetTokenT<reco::BeamSpot> FastPrimaryVertexWithWeightsProducer::beamSpotToken
private

◆ clustersToken

edm::EDGetTokenT<SiPixelClusterCollectionNew> FastPrimaryVertexWithWeightsProducer::clustersToken
private

◆ jetsToken

edm::EDGetTokenT<edm::View<reco::Jet> > FastPrimaryVertexWithWeightsProducer::jetsToken
private

◆ m_barrel

const bool FastPrimaryVertexWithWeightsProducer::m_barrel
private

Definition at line 75 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_clusters

const edm::InputTag FastPrimaryVertexWithWeightsProducer::m_clusters
private

Definition at line 65 of file FastPrimaryVertexWithWeightsProducer.cc.

◆ m_EC_weight

const double FastPrimaryVertexWithWeightsProducer::m_EC_weight
private

Definition at line 104 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_endCap

const bool FastPrimaryVertexWithWeightsProducer::m_endCap
private

Definition at line 98 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_maxDeltaPhi

const double FastPrimaryVertexWithWeightsProducer::m_maxDeltaPhi
private

Definition at line 77 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_maxDeltaPhi_EC

const double FastPrimaryVertexWithWeightsProducer::m_maxDeltaPhi_EC
private

Definition at line 101 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_maxJetEta

const double FastPrimaryVertexWithWeightsProducer::m_maxJetEta
private

Definition at line 73 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_maxJetEta_EC

const double FastPrimaryVertexWithWeightsProducer::m_maxJetEta_EC
private

Definition at line 100 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_maxSizeX

const double FastPrimaryVertexWithWeightsProducer::m_maxSizeX
private

Definition at line 76 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_maxSizeY_q

const double FastPrimaryVertexWithWeightsProducer::m_maxSizeY_q
private

Definition at line 85 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_maxZ

const double FastPrimaryVertexWithWeightsProducer::m_maxZ
private

Definition at line 64 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_minJetEta_EC

const double FastPrimaryVertexWithWeightsProducer::m_minJetEta_EC
private

Definition at line 99 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_minJetPt

const double FastPrimaryVertexWithWeightsProducer::m_minJetPt
private

Definition at line 74 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_minSizeY_q

const double FastPrimaryVertexWithWeightsProducer::m_minSizeY_q
private

Definition at line 83 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_njets

const int FastPrimaryVertexWithWeightsProducer::m_njets
private

Definition at line 72 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_peakSizeY_q

const double FastPrimaryVertexWithWeightsProducer::m_peakSizeY_q
private

Definition at line 95 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_PixelCellHeightOverWidth

const double FastPrimaryVertexWithWeightsProducer::m_PixelCellHeightOverWidth
private

Definition at line 81 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_pixelCPE

std::string FastPrimaryVertexWithWeightsProducer::m_pixelCPE
private

Definition at line 66 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_ptWeighting

const bool FastPrimaryVertexWithWeightsProducer::m_ptWeighting
private

Definition at line 122 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_ptWeighting_offset

const double FastPrimaryVertexWithWeightsProducer::m_ptWeighting_offset
private

Definition at line 124 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_ptWeighting_slope

const double FastPrimaryVertexWithWeightsProducer::m_ptWeighting_slope
private

Definition at line 123 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_weight_charge_down

const double FastPrimaryVertexWithWeightsProducer::m_weight_charge_down
private

Definition at line 78 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_weight_charge_peak

const double FastPrimaryVertexWithWeightsProducer::m_weight_charge_peak
private

Definition at line 93 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_weight_charge_up

const double FastPrimaryVertexWithWeightsProducer::m_weight_charge_up
private

Definition at line 79 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_weight_dPhi

const double FastPrimaryVertexWithWeightsProducer::m_weight_dPhi
private

Definition at line 90 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_weight_dPhi_EC

const double FastPrimaryVertexWithWeightsProducer::m_weight_dPhi_EC
private

Definition at line 105 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_weight_rho_up

const double FastPrimaryVertexWithWeightsProducer::m_weight_rho_up
private

Definition at line 92 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_weight_SizeX1

const double FastPrimaryVertexWithWeightsProducer::m_weight_SizeX1
private

Definition at line 91 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_weightCut_step2

const double FastPrimaryVertexWithWeightsProducer::m_weightCut_step2
private

Definition at line 114 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_weightCut_step3

const double FastPrimaryVertexWithWeightsProducer::m_weightCut_step3
private

Definition at line 119 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_zClusterSearchArea_step2

const double FastPrimaryVertexWithWeightsProducer::m_zClusterSearchArea_step2
private

Definition at line 113 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_zClusterSearchArea_step3

const double FastPrimaryVertexWithWeightsProducer::m_zClusterSearchArea_step3
private

Definition at line 118 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_zClusterWidth_step1

const double FastPrimaryVertexWithWeightsProducer::m_zClusterWidth_step1
private

Definition at line 109 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_zClusterWidth_step2

const double FastPrimaryVertexWithWeightsProducer::m_zClusterWidth_step2
private

Definition at line 112 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

◆ m_zClusterWidth_step3

const double FastPrimaryVertexWithWeightsProducer::m_zClusterWidth_step3
private

Definition at line 117 of file FastPrimaryVertexWithWeightsProducer.cc.

Referenced by produce().

GeomDet::position
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:43
TrackerGeometry::idToDet
const TrackerGeomDet * idToDet(DetId) const override
Definition: TrackerGeometry.cc:193
edm::ESHandle::product
T const * product() const
Definition: ESHandle.h:86
PixelClusterParameterEstimator
Definition: PixelClusterParameterEstimator.h:15
FastPrimaryVertexWithWeightsProducer::m_maxZ
const double m_maxZ
Definition: FastPrimaryVertexWithWeightsProducer.cc:64
FastPrimaryVertexWithWeightsProducer::m_zClusterWidth_step2
const double m_zClusterWidth_step2
Definition: FastPrimaryVertexWithWeightsProducer.cc:112
mps_fire.i
i
Definition: mps_fire.py:355
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
pwdgSkimBPark_cfi.beamSpot
beamSpot
Definition: pwdgSkimBPark_cfi.py:5
LumiMonitor_cff.pixelClusters
pixelClusters
Definition: LumiMonitor_cff.py:10
edm::Handle::product
T const * product() const
Definition: Handle.h:70
FastPrimaryVertexWithWeightsProducer::m_minSizeY_q
const double m_minSizeY_q
Definition: FastPrimaryVertexWithWeightsProducer.cc:83
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
multPhiCorr_741_25nsDY_cfi.py
py
Definition: multPhiCorr_741_25nsDY_cfi.py:12
edm
HLT enums.
Definition: AlignableModifier.h:19
reco::Vertex::Error
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
edmNew::DetSetVector::const_iterator
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetSetVectorNew.h:231
TkPixelCPERecord
Definition: TkPixelCPERecord.h:18
FastPrimaryVertexWithWeightsProducer::m_barrel
const bool m_barrel
Definition: FastPrimaryVertexWithWeightsProducer.cc:75
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
SiPixelCluster
Pixel cluster – collection of neighboring pixels above threshold.
Definition: SiPixelCluster.h:27
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
findQualityFiles.v
v
Definition: findQualityFiles.py:179
FastPrimaryVertexWithWeightsProducer::m_zClusterSearchArea_step2
const double m_zClusterSearchArea_step2
Definition: FastPrimaryVertexWithWeightsProducer.cc:113
edm::Handle
Definition: AssociativeIterator.h:50
FastPrimaryVertexWithWeightsProducer::m_endCap
const bool m_endCap
Definition: FastPrimaryVertexWithWeightsProducer.cc:98
FastPrimaryVertexWithWeightsProducer::clustersToken
edm::EDGetTokenT< SiPixelClusterCollectionNew > clustersToken
Definition: FastPrimaryVertexWithWeightsProducer.cc:67
FastPrimaryVertexWithWeightsProducer::m_weightCut_step3
const double m_weightCut_step3
Definition: FastPrimaryVertexWithWeightsProducer.cc:119
edmNew::DetSet::size
size_type size() const
Definition: DetSetNew.h:68
TrackerGeometry::idToDetUnit
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: TrackerGeometry.cc:183
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
DetId
Definition: DetId.h:17
GeomDet::surface
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
SiPixelCluster::sizeY
int sizeY() const
Definition: SiPixelCluster.h:128
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
FastPrimaryVertexWithWeightsProducer::m_maxDeltaPhi
const double m_maxDeltaPhi
Definition: FastPrimaryVertexWithWeightsProducer.cc:77
FastPrimaryVertexWithWeightsProducer::m_PixelCellHeightOverWidth
const double m_PixelCellHeightOverWidth
Definition: FastPrimaryVertexWithWeightsProducer.cc:81
FastPrimaryVertexWithWeightsProducer::m_maxSizeY_q
const double m_maxSizeY_q
Definition: FastPrimaryVertexWithWeightsProducer.cc:85
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
SiPixelRawToDigiRegional_cfi.deltaPhi
deltaPhi
Definition: SiPixelRawToDigiRegional_cfi.py:9
PVValHelper::eta
Definition: PVValidationHelpers.h:69
TrackerDigiGeometryRecord
Definition: TrackerDigiGeometryRecord.h:15
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
Surface::toGlobal
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
DDAxes::z
edm::ESHandle
Definition: DTSurvey.h:22
edmNew::DetSet
Definition: DetSetNew.h:22
RecoTauValidation_cfi.sizeY
sizeY
Definition: RecoTauValidation_cfi.py:291
FastPrimaryVertexWithWeightsProducer::m_zClusterSearchArea_step3
const double m_zClusterSearchArea_step3
Definition: FastPrimaryVertexWithWeightsProducer.cc:118
Point3DBase< float, GlobalTag >
FastPrimaryVertexWithWeightsProducer::m_pixelCPE
std::string m_pixelCPE
Definition: FastPrimaryVertexWithWeightsProducer.cc:66
DDAxes::rho
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
FastPrimaryVertexWithWeightsProducer::m_weight_dPhi_EC
const double m_weight_dPhi_EC
Definition: FastPrimaryVertexWithWeightsProducer.cc:105
PbPb_ZMuSkimMuonDPG_cff.tracker
tracker
Definition: PbPb_ZMuSkimMuonDPG_cff.py:60
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
edm::View
Definition: CaloClusterFwd.h:14
FastPrimaryVertexWithWeightsProducer::m_zClusterWidth_step1
const double m_zClusterWidth_step1
Definition: FastPrimaryVertexWithWeightsProducer.cc:109
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
FastPrimaryVertexWithWeightsProducer::m_minJetPt
const double m_minJetPt
Definition: FastPrimaryVertexWithWeightsProducer.cc:74
FastPrimaryVertexWithWeightsProducer::m_maxSizeX
const double m_maxSizeX
Definition: FastPrimaryVertexWithWeightsProducer.cc:76
HLT_2018_cff.weight_dPhi
weight_dPhi
Definition: HLT_2018_cff.py:49259
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
iEvent
int iEvent
Definition: GenABIO.cc:224
FastPrimaryVertexWithWeightsProducer::m_maxJetEta_EC
const double m_maxJetEta_EC
Definition: FastPrimaryVertexWithWeightsProducer.cc:100
FastPrimaryVertexWithWeightsProducer::beamSpotToken
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken
Definition: FastPrimaryVertexWithWeightsProducer.cc:68
FastPrimaryVertexWithWeightsProducer::m_ptWeighting
const bool m_ptWeighting
Definition: FastPrimaryVertexWithWeightsProducer.cc:122
get
#define get
FastPrimaryVertexWithWeightsProducer::m_weight_rho_up
const double m_weight_rho_up
Definition: FastPrimaryVertexWithWeightsProducer.cc:92
res
Definition: Electron.h:6
FastPrimaryVertexWithWeightsProducer::m_maxDeltaPhi_EC
const double m_maxDeltaPhi_EC
Definition: FastPrimaryVertexWithWeightsProducer.cc:101
FastPrimaryVertexWithWeightsProducer::m_weight_charge_down
const double m_weight_charge_down
Definition: FastPrimaryVertexWithWeightsProducer.cc:78
FastPrimaryVertexWithWeightsProducer::m_maxJetEta
const double m_maxJetEta
Definition: FastPrimaryVertexWithWeightsProducer.cc:73
FastPrimaryVertexWithWeightsProducer::jetsToken
edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken
Definition: FastPrimaryVertexWithWeightsProducer.cc:69
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
FastPrimaryVertexWithWeightsProducer::m_ptWeighting_offset
const double m_ptWeighting_offset
Definition: FastPrimaryVertexWithWeightsProducer.cc:124
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
edmNew::DetSetVector
Definition: DetSetNew.h:13
FindPeakFastPV
float FindPeakFastPV(const std::vector< float > &zProjections, const std::vector< float > &zWeights, const float oldVertex, const float m_zClusterWidth, const float m_zClusterSearchArea, const float m_weightCut)
Definition: FindPeakFastPV.h:13
reco::Vertex::Point
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
FastPrimaryVertexWithWeightsProducer::m_weight_dPhi
const double m_weight_dPhi
Definition: FastPrimaryVertexWithWeightsProducer.cc:90
FastPrimaryVertexWithWeightsProducer::m_EC_weight
const double m_EC_weight
Definition: FastPrimaryVertexWithWeightsProducer.cc:104
FastPrimaryVertexWithWeightsProducer::m_minJetEta_EC
const double m_minJetEta_EC
Definition: FastPrimaryVertexWithWeightsProducer.cc:99
SiPixelCluster::charge
int charge() const
Definition: SiPixelCluster.h:130
FastPrimaryVertexWithWeightsProducer::m_njets
const int m_njets
Definition: FastPrimaryVertexWithWeightsProducer.cc:72
edm::View::const_iterator
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
FastPrimaryVertexWithWeightsProducer::m_weight_SizeX1
const double m_weight_SizeX1
Definition: FastPrimaryVertexWithWeightsProducer.cc:91
objectSelection_cff.selectedJets
selectedJets
Definition: objectSelection_cff.py:99
FastPrimaryVertexWithWeightsProducer::m_weight_charge_up
const double m_weight_charge_up
Definition: FastPrimaryVertexWithWeightsProducer.cc:79
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
createTree.pp
pp
Definition: createTree.py:17
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
SiPixelCluster::sizeX
int sizeX() const
Definition: SiPixelCluster.h:125
FastPrimaryVertexWithWeightsProducer::m_peakSizeY_q
const double m_peakSizeY_q
Definition: FastPrimaryVertexWithWeightsProducer.cc:95
FastPrimaryVertexWithWeightsProducer::m_ptWeighting_slope
const double m_ptWeighting_slope
Definition: FastPrimaryVertexWithWeightsProducer.cc:123
edm::InputTag
Definition: InputTag.h:15
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
weight
Definition: weight.py:1
reco::Vertex
Definition: Vertex.h:35
FastPrimaryVertexWithWeightsProducer::m_weightCut_step2
const double m_weightCut_step2
Definition: FastPrimaryVertexWithWeightsProducer.cc:114
FastPrimaryVertexWithWeightsProducer::m_zClusterWidth_step3
const double m_zClusterWidth_step3
Definition: FastPrimaryVertexWithWeightsProducer.cc:117
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
TrackerGeometry
Definition: TrackerGeometry.h:14
FastPrimaryVertexWithWeightsProducer::m_weight_charge_peak
const double m_weight_charge_peak
Definition: FastPrimaryVertexWithWeightsProducer.cc:93