![]() |
![]() |
#include <RecoBTag/FastPrimaryVertexProducer/src/FastPrimaryVertexProducer.cc>
Public Member Functions | |
FastPrimaryVertexProducer (const edm::ParameterSet &) | |
Private Member Functions | |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
Private Attributes | |
edm::InputTag | m_beamSpot |
double | m_clusterLength |
edm::InputTag | m_clusters |
edm::InputTag | m_jets |
double | m_maxDeltaPhi |
double | m_maxSizeX |
double | m_maxZ |
std::string | m_pixelCPE |
Description: [one line class summary]
Implementation: [Notes on implementation]
Definition at line 75 of file FastPrimaryVertexProducer.cc.
FastPrimaryVertexProducer::FastPrimaryVertexProducer | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 91 of file FastPrimaryVertexProducer.cc.
References edm::ParameterSet::getParameter().
{ m_clusters = iConfig.getParameter<edm::InputTag>("clusters"); m_jets = iConfig.getParameter<edm::InputTag>("jets"); m_beamSpot = iConfig.getParameter<edm::InputTag>("beamSpot"); m_pixelCPE = iConfig.getParameter<std::string>("pixelCPE"); m_maxZ = iConfig.getParameter<double>("maxZ"); m_maxSizeX = iConfig.getParameter<double>("maxSizeX"); m_maxDeltaPhi = iConfig.getParameter<double>("maxDeltaPhi"); m_clusterLength = iConfig.getParameter<double>("clusterLength"); produces<reco::VertexCollection>(); }
void FastPrimaryVertexProducer::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 107 of file FastPrimaryVertexProducer.cc.
References SiPixelRawToDigiRegional_cfi::beamSpot, SiPixelRawToDigiRegional_cfi::deltaPhi, alignCSCRings::e, edm::EventSetup::get(), edm::Event::getByLabel(), j, analyzePatCleaning_cfg::jets, max(), n, AlCaHLTBitMon_ParallelJobs::p, PV3DBase< T, PVType, FrameType >::phi(), clustersummaryproducer_cfg::pixelClusters, createTree::pp, edm::Event::put(), dt_dqm_sourceclient_common_cff::reco, edmNew::DetSet< T >::size(), SiPixelCluster::sizeX(), SiPixelCluster::sizeY(), python::multivaluedict::sort(), patCandidatesForDimuonsSequences_cff::tracker, v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and z.
{ using namespace edm; using namespace reco; using namespace std; Handle<SiPixelClusterCollectionNew> cH; iEvent.getByLabel(m_clusters,cH); const SiPixelClusterCollectionNew & pixelClusters = *cH.product(); Handle<edm::View<reco::Jet> > jH; iEvent.getByLabel(m_jets,jH); const edm::View<reco::Jet> & jets = *jH.product(); CaloJetCollection selectedJets; for(edm::View<reco::Jet>::const_iterator it = jets.begin() ; it != jets.end() ; it++) { if(it->pt() > 40 && fabs(it->eta()) < 1.6) { const CaloJet * ca = dynamic_cast<const CaloJet *>( &(*it)); if(ca ==0) abort(); selectedJets.push_back(*ca); // std::cout << "Jet eta,phi,pt: "<< it->eta() << "," << it->phi() << "," << it->pt() << std::endl; } } edm::ESHandle<PixelClusterParameterEstimator> pe; const PixelClusterParameterEstimator * pp ; iSetup.get<TkPixelCPERecord>().get(m_pixelCPE , pe ); pp = pe.product(); edm::Handle<BeamSpot> beamSpot; iEvent.getByLabel(m_beamSpot,beamSpot); edm::ESHandle<TrackerGeometry> tracker; iSetup.get<TrackerDigiGeometryRecord>().get(tracker); const TrackerGeometry * trackerGeometry = tracker.product(); float lengthBmodule=6.66;//cm std::vector<float> zProjections; for(CaloJetCollection::const_iterator jit = selectedJets.begin() ; jit != selectedJets.end() ; jit++) { float px=jit->px(); float py=jit->py(); float pz=jit->pz(); float pt=jit->pt(); float jetZOverRho = jit->momentum().Z()/jit->momentum().Rho(); int minSizeY = fabs(2.*jetZOverRho)-1; int maxSizeY = fabs(2.*jetZOverRho)+2; if( fabs(jit->eta()) > 1.6) { minSizeY = 1; } for(SiPixelClusterCollectionNew::const_iterator it = pixelClusters.begin() ; it != pixelClusters.end() ; it++) //Loop on pixel modules with clusters { DetId id = it->detId(); const edmNew::DetSet<SiPixelCluster> & detset = (*it); Point3DBase<float, GlobalTag> modulepos=trackerGeometry->idToDet(id)->position(); float zmodule = modulepos.z() - ((modulepos.x()-beamSpot->x0())*px+(modulepos.y()-beamSpot->y0())*py)/pt * pz/pt; if ((fabs(deltaPhi(jit->momentum().Phi(),modulepos.phi()))< m_maxDeltaPhi*2)&&(fabs(zmodule)<(m_maxZ+lengthBmodule/2))){ for(size_t j = 0 ; j < detset.size() ; j ++) // Loop on pixel clusters on this module { const SiPixelCluster & aCluster = detset[j]; if(aCluster.sizeX() < m_maxSizeX && aCluster.sizeY() >= minSizeY && aCluster.sizeY() <= maxSizeY) { Point3DBase<float, GlobalTag> v = trackerGeometry->idToDet(id)->surface().toGlobal(pp->localParametersV( aCluster,( *trackerGeometry->idToDetUnit(id)))[0].first) ; GlobalPoint v_bs(v.x()-beamSpot->x0(),v.y()-beamSpot->y0(),v.z()); if(fabs(deltaPhi(jit->momentum().Phi(),v_bs.phi())) < m_maxDeltaPhi) { float z = v.z() - ((v.x()-beamSpot->x0())*px+(v.y()-beamSpot->y0())*py)/pt * pz/pt; if(fabs(z) < m_maxZ) { zProjections.push_back(z); } } } //if compatible cluster } // loop on module hits } // if compatible module } // loop on pixel modules } // loop on selected jets std::sort(zProjections.begin(),zProjections.end()); std::vector<float>::iterator itCenter = zProjections.begin(); std::vector<float>::iterator itLeftSide = zProjections.begin(); std::vector<float>::iterator itRightSide = zProjections.begin(); std::vector<int> counts; float zCluster = m_clusterLength/2.0; //cm int max=0; std::vector<float>::iterator left,right; for(;itCenter!=zProjections.end(); itCenter++) { while(itLeftSide != zProjections.end() && (*itCenter - *itLeftSide) > zCluster ) itLeftSide++; while(itRightSide != zProjections.end() && (*itRightSide - *itCenter) < zCluster ) itRightSide++; int n= itRightSide-itLeftSide; // std::cout << "algo :"<< *itCenter << " " << itCenter-zProjections.begin() << " dists: " << (*itCenter - *itLeftSide) << " " << (*itRightSide - *itCenter) << " count: " << n << std::endl; counts.push_back(n); if(n > max) { max=n; left=itLeftSide; } if(n >= max) { max=n; right=itRightSide; // std::cout << "algo :"<< i << " " << j << " " << *itCenter << " " << itCenter-zProjections.begin() << " dists: " << (*itCenter - *itLeftSide) << " " << (*itRightSide - *itCenter) << " count: " << n << std::endl; } } float res=0; if(zProjections.size() > 0) { res=*(left+(right-left)/2); // std::cout << "RES " << res << std::endl; Vertex::Error e; e(0, 0) = 0.0015 * 0.0015; e(1, 1) = 0.0015 * 0.0015; e(2, 2) = 1.5 * 1.5; Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res); Vertex thePV(p, e, 1, 1, 0); std::auto_ptr<reco::VertexCollection> pOut(new reco::VertexCollection()); pOut->push_back(thePV); iEvent.put(pOut); } else { // std::cout << "DUMMY " << res << std::endl; Vertex::Error e; e(0, 0) = 0.0015 * 0.0015; e(1, 1) = 0.0015 * 0.0015; e(2, 2) = 1.5 * 1.5; Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res); Vertex thePV(p, e, 0, 0, 0); std::auto_ptr<reco::VertexCollection> pOut(new reco::VertexCollection()); pOut->push_back(thePV); iEvent.put(pOut); } }
Definition at line 83 of file FastPrimaryVertexProducer.cc.
double FastPrimaryVertexProducer::m_clusterLength [private] |
Definition at line 88 of file FastPrimaryVertexProducer.cc.
Definition at line 81 of file FastPrimaryVertexProducer.cc.
Definition at line 82 of file FastPrimaryVertexProducer.cc.
double FastPrimaryVertexProducer::m_maxDeltaPhi [private] |
Definition at line 87 of file FastPrimaryVertexProducer.cc.
double FastPrimaryVertexProducer::m_maxSizeX [private] |
Definition at line 86 of file FastPrimaryVertexProducer.cc.
double FastPrimaryVertexProducer::m_maxZ [private] |
Definition at line 85 of file FastPrimaryVertexProducer.cc.
std::string FastPrimaryVertexProducer::m_pixelCPE [private] |
Definition at line 84 of file FastPrimaryVertexProducer.cc.