CMS 3D CMS Logo

FastPrimaryVertexProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FastPrimaryVertexProducer
4 // Class: FastPrimaryVertexProducer
5 //
13 //
14 // Original Author: Andrea RIZZI
15 // Created: Thu Dec 22 14:51:44 CET 2011
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 #include <vector>
22 // user include files
28 
33 
36 
46 
50 
54 #define HaveMtv
55 #define HaveFsmw
56 #define HaveDivisive
57 #ifdef HaveMtv
59 #endif
60 #ifdef HaveFsmw
62 #endif
63 #ifdef HaveDivisive
65 #endif
66 
67 using namespace std;
68 
69 //
70 // class declaration
71 //
72 
74 public:
76 
77 private:
78  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
83  double m_maxZ;
84  double m_maxSizeX;
85  double m_maxDeltaPhi;
87 };
88 
90  m_clusters = consumes<SiPixelClusterCollectionNew>(iConfig.getParameter<edm::InputTag>("clusters"));
91  m_jets = consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("jets"));
92  m_beamSpot = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"));
93  m_pixelCPE = iConfig.getParameter<std::string>("pixelCPE");
94  m_maxZ = iConfig.getParameter<double>("maxZ");
95  m_maxSizeX = iConfig.getParameter<double>("maxSizeX");
96  m_maxDeltaPhi = iConfig.getParameter<double>("maxDeltaPhi");
97  m_clusterLength = iConfig.getParameter<double>("clusterLength");
98  produces<reco::VertexCollection>();
99 }
100 
102  using namespace edm;
103  using namespace reco;
104  using namespace std;
105 
107  iEvent.getByToken(m_clusters, cH);
109 
111  iEvent.getByToken(m_jets, jH);
112  const edm::View<reco::Jet>& jets = *jH.product();
113 
115  for (edm::View<reco::Jet>::const_iterator it = jets.begin(); it != jets.end(); it++) {
116  if (it->pt() > 40 && fabs(it->eta()) < 1.6) {
117  const CaloJet* ca = dynamic_cast<const CaloJet*>(&(*it));
118  if (ca == nullptr)
119  abort();
120  selectedJets.push_back(*ca);
121  // std::cout << "Jet eta,phi,pt: "<< it->eta() << "," << it->phi() << "," << it->pt() << std::endl;
122  }
123  }
124 
127  iSetup.get<TkPixelCPERecord>().get(m_pixelCPE, pe);
128  pp = pe.product();
129 
131  iEvent.getByToken(m_beamSpot, beamSpot);
132 
134  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
135  const TrackerGeometry* trackerGeometry = tracker.product();
136 
137  float lengthBmodule = 6.66; //cm
138  std::vector<float> zProjections;
139  for (CaloJetCollection::const_iterator jit = selectedJets.begin(); jit != selectedJets.end(); jit++) {
140  float px = jit->px();
141  float py = jit->py();
142  float pz = jit->pz();
143  float pt = jit->pt();
144 
145  float jetZOverRho = jit->momentum().Z() / jit->momentum().Rho();
146  int minSizeY = fabs(2. * jetZOverRho) - 1;
147  int maxSizeY = fabs(2. * jetZOverRho) + 2;
148  if (fabs(jit->eta()) > 1.6) {
149  minSizeY = 1;
150  }
151 
152  for (SiPixelClusterCollectionNew::const_iterator it = pixelClusters.begin(); it != pixelClusters.end();
153  it++) //Loop on pixel modules with clusters
154  {
155  DetId id = it->detId();
156  const edmNew::DetSet<SiPixelCluster>& detset = (*it);
157  Point3DBase<float, GlobalTag> modulepos = trackerGeometry->idToDet(id)->position();
158  float zmodule = modulepos.z() -
159  ((modulepos.x() - beamSpot->x0()) * px + (modulepos.y() - beamSpot->y0()) * py) / pt * pz / pt;
160  if ((fabs(deltaPhi(jit->momentum().Phi(), modulepos.phi())) < m_maxDeltaPhi * 2) &&
161  (fabs(zmodule) < (m_maxZ + lengthBmodule / 2))) {
162  for (size_t j = 0; j < detset.size(); j++) // Loop on pixel clusters on this module
163  {
164  const SiPixelCluster& aCluster = detset[j];
165  if (aCluster.sizeX() < m_maxSizeX && aCluster.sizeY() >= minSizeY && aCluster.sizeY() <= maxSizeY) {
166  Point3DBase<float, GlobalTag> v = trackerGeometry->idToDet(id)->surface().toGlobal(
167  pp->localParametersV(aCluster, (*trackerGeometry->idToDetUnit(id)))[0].first);
168  GlobalPoint v_bs(v.x() - beamSpot->x0(), v.y() - beamSpot->y0(), v.z());
169  if (fabs(deltaPhi(jit->momentum().Phi(), v_bs.phi())) < m_maxDeltaPhi) {
170  float z = v.z() - ((v.x() - beamSpot->x0()) * px + (v.y() - beamSpot->y0()) * py) / pt * pz / pt;
171  if (fabs(z) < m_maxZ) {
172  zProjections.push_back(z);
173  }
174  }
175  } //if compatible cluster
176  } // loop on module hits
177  } // if compatible module
178  } // loop on pixel modules
179 
180  } // loop on selected jets
181  std::sort(zProjections.begin(), zProjections.end());
182 
183  std::vector<float>::iterator itCenter = zProjections.begin();
184  std::vector<float>::iterator itLeftSide = zProjections.begin();
185  std::vector<float>::iterator itRightSide = zProjections.begin();
186  std::vector<int> counts;
187  float zCluster = m_clusterLength / 2.0; //cm
188  int max = 0;
189  std::vector<float>::iterator left, right;
190  for (; itCenter != zProjections.end(); itCenter++) {
191  while (itLeftSide != zProjections.end() && (*itCenter - *itLeftSide) > zCluster)
192  itLeftSide++;
193  while (itRightSide != zProjections.end() && (*itRightSide - *itCenter) < zCluster)
194  itRightSide++;
195 
196  int n = itRightSide - itLeftSide;
197  // std::cout << "algo :"<< *itCenter << " " << itCenter-zProjections.begin() << " dists: " << (*itCenter - *itLeftSide) << " " << (*itRightSide - *itCenter) << " count: " << n << std::endl;
198  counts.push_back(n);
199  if (n > max) {
200  max = n;
201  left = itLeftSide;
202  }
203  if (n >= max) {
204  max = n;
205  right = itRightSide;
206  // std::cout << "algo :"<< i << " " << j << " " << *itCenter << " " << itCenter-zProjections.begin() << " dists: " << (*itCenter - *itLeftSide) << " " << (*itRightSide - *itCenter) << " count: " << n << std::endl;
207  }
208  }
209 
210  float res = 0;
211  if (!zProjections.empty()) {
212  res = *(left + (right - left) / 2);
213  // std::cout << "RES " << res << std::endl;
215  e(0, 0) = 0.0015 * 0.0015;
216  e(1, 1) = 0.0015 * 0.0015;
217  e(2, 2) = 1.5 * 1.5;
218  Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res);
219  Vertex thePV(p, e, 1, 1, 0);
220  auto pOut = std::make_unique<reco::VertexCollection>();
221  pOut->push_back(thePV);
222  iEvent.put(std::move(pOut));
223  } else {
224  // std::cout << "DUMMY " << res << std::endl;
225 
227  e(0, 0) = 0.0015 * 0.0015;
228  e(1, 1) = 0.0015 * 0.0015;
229  e(2, 2) = 1.5 * 1.5;
230  Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res);
231  Vertex thePV(p, e, 0, 0, 0);
232  auto pOut = std::make_unique<reco::VertexCollection>();
233  pOut->push_back(thePV);
234  iEvent.put(std::move(pOut));
235  }
236 }
237 
238 //define this as a plug-in
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:81
T getParameter(std::string const &) const
const_iterator end(bool update=false) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
Jets made from CaloTowers.
Definition: CaloJet.h:27
edm::EDGetTokenT< edm::View< reco::Jet > > m_jets
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
T y() const
Definition: PV3DBase.h:60
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
Definition: Electron.h:6
edm::EDGetTokenT< SiPixelClusterCollectionNew > m_clusters
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:43
const_iterator begin() const
T z() const
Definition: PV3DBase.h:61
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
FastPrimaryVertexProducer(const edm::ParameterSet &)
edm::EDGetTokenT< reco::BeamSpot > m_beamSpot
virtual VLocalValues localParametersV(const SiPixelCluster &cluster, const GeomDetUnit &gd) const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
T const * product() const
Definition: Handle.h:69
double x(const double z) const
x coordinate of the beeam spot position at a given z value (it takes into account the dxdz slope) ...
Definition: BeamSpot.h:68
int sizeY() const
Pixel cluster – collection of neighboring pixels above threshold.
fixed size matrix
HLT enums.
double y(const double z) const
y coordinate of the beeam spot position at a given z value (it takes into account the dydz slope) ...
Definition: BeamSpot.h:70
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
T get() const
Definition: EventSetup.h:73
const TrackerGeomDet * idToDet(DetId) const override
double y0() const
y coordinate
Definition: BeamSpot.h:63
const_iterator end() const
size_type size() const
Definition: DetSetNew.h:68
int sizeX() const
T x() const
Definition: PV3DBase.h:59
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511
const_iterator begin(bool update=false) const
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
double x0() const
x coordinate
Definition: BeamSpot.h:61
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override