CMS 3D CMS Logo

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

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

Inheritance diagram for FastPrimaryVertexProducer:
edm::global::EDProducer<> edm::global::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 FastPrimaryVertexProducer (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::global::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsStreamLuminosityBlocks () const final
 
bool wantsStreamRuns () const final
 
- Public Member Functions inherited from edm::global::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () 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)
 
 ~ProducerBase () noexcept(false) override
 
- 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
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
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)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

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

Private Attributes

edm::EDGetTokenT< reco::BeamSpotm_beamSpot
 
double m_clusterLength
 
edm::EDGetTokenT< SiPixelClusterCollectionNewm_clusters
 
edm::EDGetTokenT< edm::View< reco::Jet > > m_jets
 
double m_maxDeltaPhi
 
double m_maxSizeX
 
double m_maxZ
 
std::string m_pixelCPE
 

Additional Inherited Members

- Public Types inherited from edm::global::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::global::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
ProducesCollector producesCollector ()
 
- 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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 73 of file FastPrimaryVertexProducer.cc.

Constructor & Destructor Documentation

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

Definition at line 89 of file FastPrimaryVertexProducer.cc.

References edm::ParameterSet::getParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

89  {
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 }
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::View< reco::Jet > > m_jets
edm::EDGetTokenT< SiPixelClusterCollectionNew > m_clusters
edm::EDGetTokenT< reco::BeamSpot > m_beamSpot

Member Function Documentation

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

Definition at line 101 of file FastPrimaryVertexProducer.cc.

References pwdgSkimBPark_cfi::beamSpot, edm::View< T >::begin(), edmNew::DetSetVector< T >::begin(), dqmiodumpmetadata::counts, DEFINE_FWK_MODULE, SiPixelRawToDigiRegional_cfi::deltaPhi, MillePedeFileConverter_cfg::e, edm::View< T >::end(), edmNew::DetSetVector< T >::end(), edm::EventSetup::get(), edm::Event::getByToken(), TrackerGeometry::idToDet(), TrackerGeometry::idToDetUnit(), dqmiolumiharvest::j, singleTopDQM_cfi::jets, PixelClusterParameterEstimator::localParametersV(), SiStripPI::max, eostools::move(), dqmiodumpmetadata::n, 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, edm::Event::put(), multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, objectSelection_cff::selectedJets, edmNew::DetSet< T >::size(), SiPixelCluster::sizeX(), SiPixelCluster::sizeY(), GeomDet::surface(), Surface::toGlobal(), PbPb_ZMuSkimMuonDPG_cff::tracker, findQualityFiles::v, PV3DBase< T, PVType, FrameType >::x(), reco::BeamSpot::x(), reco::BeamSpot::x0(), PV3DBase< T, PVType, FrameType >::y(), reco::BeamSpot::y(), reco::BeamSpot::y0(), and PV3DBase< T, PVType, FrameType >::z().

101  {
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 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:81
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.
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
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

Member Data Documentation

edm::EDGetTokenT<reco::BeamSpot> FastPrimaryVertexProducer::m_beamSpot
private

Definition at line 81 of file FastPrimaryVertexProducer.cc.

double FastPrimaryVertexProducer::m_clusterLength
private

Definition at line 86 of file FastPrimaryVertexProducer.cc.

edm::EDGetTokenT<SiPixelClusterCollectionNew> FastPrimaryVertexProducer::m_clusters
private

Definition at line 79 of file FastPrimaryVertexProducer.cc.

edm::EDGetTokenT<edm::View<reco::Jet> > FastPrimaryVertexProducer::m_jets
private

Definition at line 80 of file FastPrimaryVertexProducer.cc.

double FastPrimaryVertexProducer::m_maxDeltaPhi
private

Definition at line 85 of file FastPrimaryVertexProducer.cc.

double FastPrimaryVertexProducer::m_maxSizeX
private

Definition at line 84 of file FastPrimaryVertexProducer.cc.

double FastPrimaryVertexProducer::m_maxZ
private

Definition at line 83 of file FastPrimaryVertexProducer.cc.

std::string FastPrimaryVertexProducer::m_pixelCPE
private

Definition at line 82 of file FastPrimaryVertexProducer.cc.