CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastPrimaryVertexWithWeightsProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FastPrimaryVertexWithWeightsProducer
4 // Class: FastPrimaryVertexWithWeightsProducer
5 //
14 //
15 // Original Author: Silvio DONATO
16 // Created: Wed Dec 18 10:05:40 CET 2013
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 #include <vector>
23 #include <math.h>
24 // user include files
33 
44 
49 
51 
53 
54 using namespace std;
55 
57  public:
59  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
60 
61  private:
62  virtual void produce(edm::Event&, const edm::EventSetup&);
63 
64 
65 
66  const double m_maxZ; // Use only pixel clusters with |z| < maxZ
67  const edm::InputTag m_clusters; // PixelClusters InputTag
68  std::string m_pixelCPE; // PixelCPE (PixelClusterParameterEstimator)
72 
73 // PARAMETERS USED IN THE BARREL PIXEL CLUSTERS PROJECTION
74  const int m_njets; // Use only the first njets
75  const double m_maxJetEta; // Use only jets with |eta| < maxJetEta
76  const double m_minJetPt; // Use only jets with Pt > minJetPt
77  const bool m_barrel; // Use clusters from pixel endcap
78  const double m_maxSizeX; // Use only pixel clusters with sizeX <= maxSizeX
79  const double m_maxDeltaPhi; // Use only pixel clusters with DeltaPhi(Jet,Cluster) < maxDeltaPhi
80  const double m_weight_charge_down; // Use only pixel clusters with ClusterCharge > weight_charge_down
81  const double m_weight_charge_up; // Use only pixel clusters with ClusterCharge < weight_charge_up
82  const double m_PixelCellHeightOverWidth;//It is the ratio between pixel cell height and width along z coordinate about 285µm/150µm=1.9
83  const double m_minSizeY_q; // Use only pixel clusters with sizeY > PixelCellHeightOverWidth * |jetZOverRho| + minSizeY_q
84  const double m_maxSizeY_q; // Use only pixel clusters with sizeY < PixelCellHeightOverWidth * |jetZOverRho| + maxSizeY_q
85 
86 // PARAMETERS USED TO WEIGHT THE BARREL PIXEL CLUSTERS
87  // The cluster weight is defined as weight = weight_dPhi * weight_sizeY * weight_rho * weight_sizeX1 * weight_charge
88 
89  const double m_weight_dPhi; // used in weight_dPhi = exp(-|DeltaPhi(JetCluster)|/m_weight_dPhi)
90  const double m_weight_SizeX1; // used in weight_SizeX1 = (ClusterSizeX==2)*1+(ClusterSizeX==1)*m_weight_SizeX1;
91  const double m_weight_rho_up; // used in weight_rho = (m_weight_rho_up - ClusterRho)/m_weight_rho_up
92  const double m_weight_charge_peak; // Give the maximum weight_charge for a cluster with Charge = m_weight_charge_peak
93  const double m_peakSizeY_q; // Give the maximum weight_sizeY for a cluster with sizeY = PixelCellHeightOverWidth * |jetZOverRho| + peakSizeY_q
94 
95 // PARAMETERS USED IN THE ENDCAP PIXEL CLUSTERS PROJECTION
96  const bool m_endCap; // Use clusters from pixel endcap
97  const double m_minJetEta_EC; // Use only jets with |eta| > minJetEta_EC
98  const double m_maxJetEta_EC; // Use only jets with |eta| < maxJetEta_EC
99  const double m_maxDeltaPhi_EC; // Use only pixel clusters with DeltaPhi(Jet,Cluster) < maxDeltaPhi_EC
100 
101 // PARAMETERS USED TO WEIGHT THE ENDCAP PIXEL CLUSTERS
102  const double m_EC_weight; // In EndCap the weight is defined as weight = m_EC_weight*(weight_dPhi)
103  const double m_weight_dPhi_EC; // Used in weight_dPhi = exp(-|DeltaPhi|/m_weight_dPhi_EC )
104 
105 // PARAMETERS USED TO FIND THE FASTPV AS PEAK IN THE Z-PROJECTIONS DISTRIBUTION
106  // First Iteration: look for a cluster with a width = m_zClusterWidth_step1
107  const double m_zClusterWidth_step1; // cluster width in step1
108 
109  // Second Iteration: use only z-projections with weight > weightCut_step2 and look for a cluster with a width = m_zClusterWidth_step2, within of weightCut_step2 of the previous result
110  const double m_zClusterWidth_step2; // cluster width in step2
111  const double m_zClusterSearchArea_step2; // cluster width in step2
112  const double m_weightCut_step2; // minimum z-projections weight required in step2
113 
114  // Third Iteration: use only z-projections with weight > weightCut_step3 and look for a cluster with a width = m_zClusterWidth_step3, within of weightCut_step3 of the previous result
115  const double m_zClusterWidth_step3; // cluster width in step3
116  const double m_zClusterSearchArea_step3; // cluster width in step3
117  const double m_weightCut_step3; // minimum z-projections weight required in step3
118 
119  // use the jetPt weighting
120  const bool m_ptWeighting; // use weight=weight*pt_weigth ?;
121  const double m_ptWeighting_slope; // pt_weigth= pt*ptWeighting_slope +m_ptWeighting_offset;
122  const double m_ptWeighting_offset; //
123 };
124 
126  m_maxZ(iConfig.getParameter<double>("maxZ")),
127  m_pixelCPE(iConfig.getParameter<std::string>("pixelCPE")),
128 
129  m_njets(iConfig.getParameter<int>("njets")),
130  m_maxJetEta(iConfig.getParameter<double>("maxJetEta")),
131  m_minJetPt(iConfig.getParameter<double>("minJetPt")),
132 
133  m_barrel(iConfig.getParameter<bool>("barrel")),
134  m_maxSizeX(iConfig.getParameter<double>("maxSizeX")),
135  m_maxDeltaPhi(iConfig.getParameter<double>("maxDeltaPhi")),
136  m_weight_charge_down(iConfig.getParameter<double>("weight_charge_down")),
137  m_weight_charge_up(iConfig.getParameter<double>("weight_charge_up")),
138  m_PixelCellHeightOverWidth(iConfig.getParameter<double>("PixelCellHeightOverWidth")),
139  m_minSizeY_q(iConfig.getParameter<double>("minSizeY_q")),
140  m_maxSizeY_q(iConfig.getParameter<double>("maxSizeY_q")),
141 
142  m_weight_dPhi(iConfig.getParameter<double>("weight_dPhi")),
143  m_weight_SizeX1(iConfig.getParameter<double>("weight_SizeX1")),
144  m_weight_rho_up(iConfig.getParameter<double>("weight_rho_up")),
145  m_weight_charge_peak(iConfig.getParameter<double>("weight_charge_peak")),
146  m_peakSizeY_q(iConfig.getParameter<double>("peakSizeY_q")),
147 
148  m_endCap(iConfig.getParameter<bool>("endCap")),
149  m_minJetEta_EC(iConfig.getParameter<double>("minJetEta_EC")),
150  m_maxJetEta_EC(iConfig.getParameter<double>("maxJetEta_EC")),
151  m_maxDeltaPhi_EC(iConfig.getParameter<double>("maxDeltaPhi_EC")),
152  m_EC_weight(iConfig.getParameter<double>("EC_weight")),
153  m_weight_dPhi_EC(iConfig.getParameter<double>("weight_dPhi_EC")),
154 
155  m_zClusterWidth_step1(iConfig.getParameter<double>("zClusterWidth_step1")),
156 
157  m_zClusterWidth_step2(iConfig.getParameter<double>("zClusterWidth_step2")),
158  m_zClusterSearchArea_step2(iConfig.getParameter<double>("zClusterSearchArea_step2")),
159  m_weightCut_step2(iConfig.getParameter<double>("weightCut_step2")),
160 
161  m_zClusterWidth_step3(iConfig.getParameter<double>("zClusterWidth_step3")),
162  m_zClusterSearchArea_step3(iConfig.getParameter<double>("zClusterSearchArea_step3")),
163  m_weightCut_step3(iConfig.getParameter<double>("weightCut_step3")),
164 
165  m_ptWeighting(iConfig.getParameter<bool>("ptWeighting")),
166  m_ptWeighting_slope(iConfig.getParameter<double>("ptWeighting_slope")),
167  m_ptWeighting_offset(iConfig.getParameter<double>("ptWeighting_offset"))
168 {
169 
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 
177 }
178 
179 void
181 {
183  desc.add<edm::InputTag> ("clusters",edm::InputTag("hltSiPixelClusters"));
184  desc.add<edm::InputTag> ("beamSpot",edm::InputTag("hltOnlineBeamSpot"));
185  desc.add<edm::InputTag> ("jets",edm::InputTag("hltCaloJetL1FastJetCorrected"));
186  desc.add<std::string>("pixelCPE","hltESPPixelCPEGeneric");
187  desc.add<double>("maxZ",19.0);
188  desc.add<int>("njets",999);
189  desc.add<double>("maxJetEta",2.6);
190  desc.add<double>("minJetPt",40.);
191  desc.add<bool>("barrel",true);
192  desc.add<double>("maxSizeX",2.1);
193  desc.add<double>("maxDeltaPhi",0.21);
194  desc.add<double>("PixelCellHeightOverWidth",1.8);
195  desc.add<double>("weight_charge_down",11.*1000.);
196  desc.add<double>("weight_charge_up",190.*1000.);
197  desc.add<double>("maxSizeY_q",2.0);
198  desc.add<double>("minSizeY_q",-0.6);
199  desc.add<double>("weight_dPhi",0.13888888);
200  desc.add<double>("weight_SizeX1",0.88);
201  desc.add<double>("weight_rho_up",22.);
202  desc.add<double>("weight_charge_peak",22.*1000.);
203  desc.add<double>("peakSizeY_q",1.0);
204  desc.add<bool>("endCap",true);
205  desc.add<double>("minJetEta_EC",1.3);
206  desc.add<double>("maxJetEta_EC",2.6);
207  desc.add<double>("maxDeltaPhi_EC",0.14);
208  desc.add<double>("EC_weight",0.008);
209  desc.add<double>("weight_dPhi_EC",0.064516129);
210  desc.add<double>("zClusterWidth_step1",2.0);
211  desc.add<double>("zClusterWidth_step2",0.65);
212  desc.add<double>("zClusterSearchArea_step2",3.0);
213  desc.add<double>("weightCut_step2",0.05);
214  desc.add<double>("zClusterWidth_step3",0.3);
215  desc.add<double>("zClusterSearchArea_step3",0.55);
216  desc.add<double>("weightCut_step3",0.1);
217  desc.add<bool>("ptWeighting",false); // <---- newMethod
218  desc.add<double>("ptWeighting_slope",1/20.);
219  desc.add<double>("ptWeighting_offset",-1);
220  descriptions.add("fastPrimaryVertexWithWeightsProducer",desc);
221 }
222 
223 void
225 {
226 
227  using namespace edm;
228  using namespace reco;
229  using namespace std;
230 
231  const float barrel_lenght=30; //half lenght of the pixel barrel 30 cm
232 
233  //get pixel cluster
235  iEvent.getByToken(clustersToken,cH);
236  const SiPixelClusterCollectionNew & pixelClusters = *cH.product();
237 
238  //get jets
240  iEvent.getByToken(jetsToken,jH);
241  const edm::View<reco::Jet> & jets = *jH.product();
242 
243  vector<const reco::Jet*> selectedJets;
244  int countjet=0;
245  for(edm::View<reco::Jet>::const_iterator it = jets.begin() ; it != jets.end() && countjet<m_njets ; it++)
246  {
247  if( //select jets used in barrel or endcap pixel cluster projections
248  ((it->pt() >= m_minJetPt) && std::abs(it->eta()) <= m_maxJetEta) || //barrel
249  ((it->pt() >= m_minJetPt) && std::abs(it->eta()) <= m_maxJetEta_EC && std::abs(it->eta()) >= m_minJetEta_EC) //endcap
250  )
251  {
252  selectedJets.push_back(&(*it));
253  countjet++;
254  }
255  }
256 
257  //get PixelClusterParameterEstimator
260  iSetup.get<TkPixelCPERecord>().get(m_pixelCPE , pe );
261  pp = pe.product();
262 
263  //get beamSpot
266 
267  //get TrackerGeometry
269  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
270  const TrackerGeometry * trackerGeometry = tracker.product();
271 
272 // PART I: get z-projections with z-weights
273  std::vector<float> zProjections;
274  std::vector<float> zWeights;
275  int jet_count=0;
276  for(vector<const reco::Jet*>::iterator jit = selectedJets.begin() ; jit != selectedJets.end() ; jit++)
277  {//loop on selected jets
278  float px=(*jit)->px();
279  float py=(*jit)->py();
280  float pz=(*jit)->pz();
281  float pt=(*jit)->pt();
282  float eta=(*jit)->eta();
283  float jetZOverRho = (*jit)->momentum().Z()/(*jit)->momentum().Rho();
284  float pt_weight= pt*m_ptWeighting_slope +m_ptWeighting_offset;
285  for(SiPixelClusterCollectionNew::const_iterator it = pixelClusters.begin() ; it != pixelClusters.end() ; it++) //Loop on pixel modules with clusters
286  {//loop on pixel modules
287  DetId id = it->detId();
288  const edmNew::DetSet<SiPixelCluster> & detset = (*it);
289  Point3DBase<float, GlobalTag> modulepos=trackerGeometry->idToDet(id)->position();
290  float zmodule = modulepos.z() - ((modulepos.x()-beamSpot->x0())*px+(modulepos.y()-beamSpot->y0())*py)/pt * pz/pt;
291  if ((std::abs(deltaPhi((*jit)->momentum().Phi(),modulepos.phi()))< m_maxDeltaPhi*2)&&(std::abs(zmodule)<(m_maxZ+barrel_lenght))){//if it is a compatible module
292  for(size_t j = 0 ; j < detset.size() ; j ++)
293  {//loop on pixel clusters on this module
294  const SiPixelCluster & aCluster = detset[j];
295  if(//it is a cluster to project
296  (// barrel
297  m_barrel &&
298  std::abs(modulepos.z())<barrel_lenght &&
299  pt>=m_minJetPt &&
300  jet_count<m_njets &&
301  std::abs(eta)<=m_maxJetEta &&
302  aCluster.sizeX() <= m_maxSizeX &&
303  aCluster.sizeY() >= m_PixelCellHeightOverWidth*std::abs(jetZOverRho)+m_minSizeY_q &&
304  aCluster.sizeY() <= m_PixelCellHeightOverWidth*std::abs(jetZOverRho)+m_maxSizeY_q
305  )
306  ||
307  (// EC
308  m_endCap &&
309  std::abs(modulepos.z())>barrel_lenght &&
310  pt > m_minJetPt &&
311  jet_count<m_njets &&
312  std::abs(eta) <= m_maxJetEta_EC &&
313  std::abs(eta) >= m_minJetEta_EC &&
314  aCluster.sizeX() <= m_maxSizeX
315  )
316  ){
317  Point3DBase<float, GlobalTag> v = trackerGeometry->idToDet(id)->surface().toGlobal(pp->localParametersV( aCluster,( *trackerGeometry->idToDetUnit(id)))[0].first) ;
318  GlobalPoint v_bs(v.x()-beamSpot->x0(),v.y()-beamSpot->y0(),v.z());
319  if(//it pass DeltaPhi(Jet,Cluster) requirements
320  (m_barrel && std::abs(modulepos.z())<barrel_lenght && std::abs(deltaPhi((*jit)->momentum().Phi(),v_bs.phi())) <= m_maxDeltaPhi ) || //barrel
321  (m_endCap && std::abs(modulepos.z())>barrel_lenght && std::abs(deltaPhi((*jit)->momentum().Phi(),v_bs.phi())) <= m_maxDeltaPhi_EC ) //EC
322  )
323  {
324  float z = v.z() - ((v.x()-beamSpot->x0())*px+(v.y()-beamSpot->y0())*py)/pt * pz/pt; //calculate z-projection
325  if(std::abs(z) < m_maxZ)
326  {
327  zProjections.push_back(z); //add z-projection in zProjections
328  float weight=0;
329  //calculate zWeight
330  if(std::abs(modulepos.z())<barrel_lenght)
331  { //barrel
332  //calculate weight_sizeY
333  float sizeY=aCluster.sizeY();
334  float sizeY_up = m_PixelCellHeightOverWidth*std::abs(jetZOverRho)+m_maxSizeY_q;
335  float sizeY_peak = m_PixelCellHeightOverWidth*std::abs(jetZOverRho)+m_peakSizeY_q;
336  float sizeY_down = m_PixelCellHeightOverWidth*std::abs(jetZOverRho)+m_minSizeY_q;
337  float weight_sizeY_up = (sizeY_up-sizeY)/(sizeY_up-sizeY_peak);
338  float weight_sizeY_down = (sizeY-sizeY_down)/(sizeY_peak-sizeY_down);
339  weight_sizeY_down = weight_sizeY_down *(weight_sizeY_down>0)*(weight_sizeY_down<1);
340  weight_sizeY_up = weight_sizeY_up *(weight_sizeY_up>0)*(weight_sizeY_up<1);
341  float weight_sizeY = weight_sizeY_up + weight_sizeY_down;
342 
343  //calculate weight_rho
344  float rho = sqrt(v_bs.x()*v_bs.x() + v_bs.y()*v_bs.y());
345  float weight_rho = ((m_weight_rho_up - rho)/m_weight_rho_up);
346 
347  //calculate weight_dPhi
348  float weight_dPhi = exp(- std::abs(deltaPhi((*jit)->momentum().Phi(),v_bs.phi()))/m_weight_dPhi );
349 
350  //calculate weight_sizeX1
351  float weight_sizeX1= (aCluster.sizeX()==2) + (aCluster.sizeX()==1)*m_weight_SizeX1;
352 
353  //calculate weight_charge
354  float charge=aCluster.charge();
356  float weightCluster_down = (charge-m_weight_charge_down)/(m_weight_charge_peak-m_weight_charge_down);
357  weightCluster_down = weightCluster_down *(weightCluster_down>0)*(weightCluster_down<1);
358  weightCluster_up = weightCluster_up *(weightCluster_up>0)*(weightCluster_up<1);
359  float weight_charge = weightCluster_up + weightCluster_down;
360 
361  //calculate the final weight
362  weight = weight_dPhi * weight_sizeY * weight_rho * weight_sizeX1 * weight_charge ;
363  }
364  else if(std::abs(modulepos.z())>barrel_lenght) // EC
365  {// EC
366  //calculate weight_dPhi
367  float weight_dPhi = exp(- std::abs(deltaPhi((*jit)->momentum().Phi(),v_bs.phi())) /m_weight_dPhi_EC );
368  //calculate the final weight
369  weight= m_EC_weight*(weight_dPhi) ;
370  }
371  if(m_ptWeighting) weight=weight*pt_weight;
372  zWeights.push_back(weight); //add the weight to zWeights
373  }
374  }//if it pass DeltaPhi(Jet,Cluster) requirements
375  }
376  }//loop on pixel clusters on this module
377  }//if it is a compatible module
378  }//loop on pixel modules
379  jet_count++;
380  }//loop on selected jets
381 
382  //order zProjections and zWeights by z
383  std::multimap<float,float> zWithW;
384  size_t i=0;
385  for(i=0;i<zProjections.size();i++) zWithW.insert(std::pair<float,float>(zProjections[i],zWeights[i]));
386  i=0;
387  for(std::multimap<float,float>::iterator it=zWithW.begin(); it!=zWithW.end(); it++,i++) { zProjections[i]=it->first; zWeights[i]=it->second; } //order zProjections and zWeights by z
388 
389 
390  //calculate zWeightsSquared
391  std::vector<float> zWeightsSquared;
392  for(std::vector<float>::iterator it=zWeights.begin();it!=zWeights.end();it++) {zWeightsSquared.push_back((*it)*(*it));}
393 
394  //do multi-step peak searching
395  float res_step1 = FindPeakFastPV( zProjections, zWeights, 0.0, m_zClusterWidth_step1, 999.0, -1.0);
396  float res_step2 = FindPeakFastPV( zProjections, zWeights, res_step1, m_zClusterWidth_step2, m_zClusterSearchArea_step2, m_weightCut_step2);
397  float res_step3 = FindPeakFastPV( zProjections, zWeightsSquared, res_step2, m_zClusterWidth_step3, m_zClusterSearchArea_step3, m_weightCut_step3*m_weightCut_step3);
398 
399  float centerWMax=res_step3;
400  //End of PART II
401 
402  //Make the output
403  float res=0;
404  if(zProjections.size() > 2)
405  {
406  res=centerWMax;
407  Vertex::Error e;
408  e(0, 0) = 0.0015 * 0.0015;
409  e(1, 1) = 0.0015 * 0.0015;
410  e(2, 2) = 1.5 * 1.5;
411  Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res);
412  Vertex thePV(p, e, 1, 1, 0);
413  std::auto_ptr<reco::VertexCollection> pOut(new reco::VertexCollection());
414  pOut->push_back(thePV);
415  iEvent.put(pOut);
416  } else
417  {
419  e(0, 0) = 0.0015 * 0.0015;
420  e(1, 1) = 0.0015 * 0.0015;
421  e(2, 2) = 1.5 * 1.5;
422  Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res);
423  Vertex thePV(p, e, 0, 0, 0);
424  std::auto_ptr<reco::VertexCollection> pOut(new reco::VertexCollection());
425  pOut->push_back(thePV);
426  iEvent.put(pOut);
427  }
428 
429 //Finally, calculate the zClusterQuality as Sum(weights near the fastPV)/sqrt(Sum(total weights)) [a kind of zCluster significance]
430 
431 
432 
433  const float half_width_peak=1;
434  float nWeightedTot=0;
435  float nWeightedTotPeak=0;
436  for(std::vector<float>::iterator it = zProjections.begin();it!=zProjections.end(); it++)
437  {
438  nWeightedTot+=zWeights[it-zProjections.begin()];
439  if((res-half_width_peak)<=(*it) && (*it)<=(res+half_width_peak))
440  {
441  nWeightedTotPeak+=zWeights[it-zProjections.begin()];
442  }
443  }
444 
445 std::auto_ptr<float > zClusterQuality(new float());
446 *zClusterQuality=-1;
447 if(nWeightedTot!=0)
448 {
449  *zClusterQuality=nWeightedTotPeak / sqrt(nWeightedTot/(2*half_width_peak)); // where 30 is the beam spot lenght
450  iEvent.put(zClusterQuality);
451 }
452 else
453  iEvent.put(zClusterQuality);
454 
455 }
456 
457 
458 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
float charge() const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
tuple pp
Definition: createTree.py:15
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T y() const
Definition: PV3DBase.h:63
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
virtual void produce(edm::Event &, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:56
void add(std::string const &label, ParameterSetDescription const &psetDescription)
tuple pixelClusters
int sizeY() const
Pixel cluster – collection of neighboring pixels above threshold.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
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)
size_type size() const
Definition: DetSetNew.h:87
int sizeX() const
edm::EDGetTokenT< SiPixelClusterCollectionNew > clustersToken
T x() const
Definition: PV3DBase.h:62