CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
 
- Public Member Functions inherited from edm::global::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducerBase ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

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

Private Attributes

edm::EDGetTokenT< reco::BeamSpotm_beamSpot
 
double m_clusterLength
 
edm::EDGetTokenT
< SiPixelClusterCollectionNew
m_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
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::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

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 74 of file FastPrimaryVertexProducer.cc.

Constructor & Destructor Documentation

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

Definition at line 91 of file FastPrimaryVertexProducer.cc.

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

92 {
93  m_clusters = consumes<SiPixelClusterCollectionNew>(iConfig.getParameter<edm::InputTag>("clusters"));
94  m_jets = consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("jets"));
95  m_beamSpot = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"));
96  m_pixelCPE = iConfig.getParameter<std::string>("pixelCPE");
97  m_maxZ = iConfig.getParameter<double>("maxZ");
98  m_maxSizeX = iConfig.getParameter<double>("maxSizeX");
99  m_maxDeltaPhi = iConfig.getParameter<double>("maxDeltaPhi");
100  m_clusterLength = iConfig.getParameter<double>("clusterLength");
101  produces<reco::VertexCollection>();
102 }
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
overrideprivatevirtual

Implements edm::global::EDProducerBase.

Definition at line 107 of file FastPrimaryVertexProducer.cc.

References SiPixelRawToDigiRegional_cfi::beamSpot, SiPixelRawToDigiRegional_cfi::deltaPhi, alignCSCRings::e, edm::EventSetup::get(), edm::Event::getByToken(), j, fwrapper::jets, bookConverter::max, gen::n, AlCaHLTBitMon_ParallelJobs::p, PV3DBase< T, PVType, FrameType >::phi(), HLT_25ns14e33_v1_cff::pixelClusters, createTree::pp, EnergyCorrector::pt, edm::Event::put(), dt_dqm_sourceclient_common_cff::reco, patRefSel_refMuJets_cfi::selectedJets, edmNew::DetSet< T >::size(), SiPixelCluster::sizeX(), SiPixelCluster::sizeY(), python.multivaluedict::sort(), patCandidatesForDimuonsSequences_cff::tracker, findQualityFiles::v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

108 {
109  using namespace edm;
110  using namespace reco;
111  using namespace std;
112 
114  iEvent.getByToken(m_clusters,cH);
115  const SiPixelClusterCollectionNew & pixelClusters = *cH.product();
116 
118  iEvent.getByToken(m_jets,jH);
119  const edm::View<reco::Jet> & jets = *jH.product();
120 
122  for(edm::View<reco::Jet>::const_iterator it = jets.begin() ; it != jets.end() ; it++)
123  {
124  if(it->pt() > 40 && fabs(it->eta()) < 1.6)
125  {
126  const CaloJet * ca = dynamic_cast<const CaloJet *>( &(*it));
127  if(ca ==0) abort();
128  selectedJets.push_back(*ca);
129 // std::cout << "Jet eta,phi,pt: "<< it->eta() << "," << it->phi() << "," << it->pt() << std::endl;
130  }
131  }
132 
135  iSetup.get<TkPixelCPERecord>().get(m_pixelCPE , pe );
136  pp = pe.product();
137 
140 
142  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
143  const TrackerGeometry * trackerGeometry = tracker.product();
144 
145 
146  float lengthBmodule=6.66;//cm
147  std::vector<float> zProjections;
148  for(CaloJetCollection::const_iterator jit = selectedJets.begin() ; jit != selectedJets.end() ; jit++)
149  {
150  float px=jit->px();
151  float py=jit->py();
152  float pz=jit->pz();
153  float pt=jit->pt();
154 
155  float jetZOverRho = jit->momentum().Z()/jit->momentum().Rho();
156  int minSizeY = fabs(2.*jetZOverRho)-1;
157  int maxSizeY = fabs(2.*jetZOverRho)+2;
158  if( fabs(jit->eta()) > 1.6)
159  {
160  minSizeY = 1;
161  }
162 
163  for(SiPixelClusterCollectionNew::const_iterator it = pixelClusters.begin() ; it != pixelClusters.end() ; it++) //Loop on pixel modules with clusters
164  {
165  DetId id = it->detId();
166  const edmNew::DetSet<SiPixelCluster> & detset = (*it);
167  Point3DBase<float, GlobalTag> modulepos=trackerGeometry->idToDet(id)->position();
168  float zmodule = modulepos.z() - ((modulepos.x()-beamSpot->x0())*px+(modulepos.y()-beamSpot->y0())*py)/pt * pz/pt;
169  if ((fabs(deltaPhi(jit->momentum().Phi(),modulepos.phi()))< m_maxDeltaPhi*2)&&(fabs(zmodule)<(m_maxZ+lengthBmodule/2))){
170 
171  for(size_t j = 0 ; j < detset.size() ; j ++) // Loop on pixel clusters on this module
172  {
173  const SiPixelCluster & aCluster = detset[j];
174  if(aCluster.sizeX() < m_maxSizeX && aCluster.sizeY() >= minSizeY && aCluster.sizeY() <= maxSizeY) {
175  Point3DBase<float, GlobalTag> v = trackerGeometry->idToDet(id)->surface().toGlobal(pp->localParametersV( aCluster,( *trackerGeometry->idToDetUnit(id)))[0].first) ;
176  GlobalPoint v_bs(v.x()-beamSpot->x0(),v.y()-beamSpot->y0(),v.z());
177  if(fabs(deltaPhi(jit->momentum().Phi(),v_bs.phi())) < m_maxDeltaPhi)
178  {
179  float z = v.z() - ((v.x()-beamSpot->x0())*px+(v.y()-beamSpot->y0())*py)/pt * pz/pt;
180  if(fabs(z) < m_maxZ)
181  {
182  zProjections.push_back(z);
183  }
184  }
185  } //if compatible cluster
186  } // loop on module hits
187  } // if compatible module
188  } // loop on pixel modules
189 
190  } // loop on selected jets
191  std::sort(zProjections.begin(),zProjections.end());
192 
193  std::vector<float>::iterator itCenter = zProjections.begin();
194  std::vector<float>::iterator itLeftSide = zProjections.begin();
195  std::vector<float>::iterator itRightSide = zProjections.begin();
196  std::vector<int> counts;
197  float zCluster = m_clusterLength/2.0; //cm
198  int max=0;
199  std::vector<float>::iterator left,right;
200  for(;itCenter!=zProjections.end(); itCenter++)
201  {
202 
203  while(itLeftSide != zProjections.end() && (*itCenter - *itLeftSide) > zCluster ) itLeftSide++;
204  while(itRightSide != zProjections.end() && (*itRightSide - *itCenter) < zCluster ) itRightSide++;
205 
206  int n= itRightSide-itLeftSide;
207  // std::cout << "algo :"<< *itCenter << " " << itCenter-zProjections.begin() << " dists: " << (*itCenter - *itLeftSide) << " " << (*itRightSide - *itCenter) << " count: " << n << std::endl;
208  counts.push_back(n);
209  if(n > max) {
210  max=n;
211  left=itLeftSide;
212  }
213  if(n >= max) {
214  max=n;
215  right=itRightSide;
216 // std::cout << "algo :"<< i << " " << j << " " << *itCenter << " " << itCenter-zProjections.begin() << " dists: " << (*itCenter - *itLeftSide) << " " << (*itRightSide - *itCenter) << " count: " << n << std::endl;
217  }
218  }
219 
220 
221 
222 
223  float res=0;
224  if(zProjections.size() > 0)
225  {
226  res=*(left+(right-left)/2);
227 // std::cout << "RES " << res << std::endl;
229  e(0, 0) = 0.0015 * 0.0015;
230  e(1, 1) = 0.0015 * 0.0015;
231  e(2, 2) = 1.5 * 1.5;
232  Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res);
233  Vertex thePV(p, e, 1, 1, 0);
234  std::auto_ptr<reco::VertexCollection> pOut(new reco::VertexCollection());
235  pOut->push_back(thePV);
236  iEvent.put(pOut);
237  } else
238  {
239  // std::cout << "DUMMY " << res << std::endl;
240 
242  e(0, 0) = 0.0015 * 0.0015;
243  e(1, 1) = 0.0015 * 0.0015;
244  e(2, 2) = 1.5 * 1.5;
245  Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res);
246  Vertex thePV(p, e, 0, 0, 0);
247  std::auto_ptr<reco::VertexCollection> pOut(new reco::VertexCollection());
248  pOut->push_back(thePV);
249  iEvent.put(pOut);
250 
251  }
252 
253 
254 }
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Jets made from CaloTowers.
Definition: CaloJet.h:29
tuple pp
Definition: createTree.py:15
edm::EDGetTokenT< edm::View< reco::Jet > > m_jets
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
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
edm::EDGetTokenT< SiPixelClusterCollectionNew > m_clusters
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
vector< PseudoJet > jets
T z() const
Definition: PV3DBase.h:64
int j
Definition: DBlmapReader.cc:9
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
edm::EDGetTokenT< reco::BeamSpot > m_beamSpot
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:56
int sizeY() const
Pixel cluster – collection of neighboring pixels above threshold.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
size_type size() const
Definition: DetSetNew.h:86
int sizeX() const
T x() const
Definition: PV3DBase.h:62
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects

Member Data Documentation

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

Definition at line 82 of file FastPrimaryVertexProducer.cc.

double FastPrimaryVertexProducer::m_clusterLength
private

Definition at line 87 of file FastPrimaryVertexProducer.cc.

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

Definition at line 80 of file FastPrimaryVertexProducer.cc.

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

Definition at line 81 of file FastPrimaryVertexProducer.cc.

double FastPrimaryVertexProducer::m_maxDeltaPhi
private

Definition at line 86 of file FastPrimaryVertexProducer.cc.

double FastPrimaryVertexProducer::m_maxSizeX
private

Definition at line 85 of file FastPrimaryVertexProducer.cc.

double FastPrimaryVertexProducer::m_maxZ
private

Definition at line 84 of file FastPrimaryVertexProducer.cc.

std::string FastPrimaryVertexProducer::m_pixelCPE
private

Definition at line 83 of file FastPrimaryVertexProducer.cc.