CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
20 // system include files
21 #include <memory>
22 #include <vector>
23 // user include files
29 
34 
37 
47 
51 
55 #define HaveMtv
56 #define HaveFsmw
57 #define HaveDivisive
58 #ifdef HaveMtv
60 #endif
61 #ifdef HaveFsmw
63 #endif
64 #ifdef HaveDivisive
66 #endif
67 
68 using namespace std;
69 
70 //
71 // class declaration
72 //
73 
75  public:
77 
78  private:
79  virtual void produce(edm::Event&, const edm::EventSetup&) override;
84  double m_maxZ;
85  double m_maxSizeX;
86  double m_maxDeltaPhi;
87  double m_clusterLength;
88 };
89 
91 {
92  m_clusters = iConfig.getParameter<edm::InputTag>("clusters");
93  m_jets = iConfig.getParameter<edm::InputTag>("jets");
94  m_beamSpot = iConfig.getParameter<edm::InputTag>("beamSpot");
95  m_pixelCPE = iConfig.getParameter<std::string>("pixelCPE");
96  m_maxZ = iConfig.getParameter<double>("maxZ");
97  m_maxSizeX = iConfig.getParameter<double>("maxSizeX");
98  m_maxDeltaPhi = iConfig.getParameter<double>("maxDeltaPhi");
99  m_clusterLength = iConfig.getParameter<double>("clusterLength");
100  produces<reco::VertexCollection>();
101 }
102 
103 
104 
105 void
107 {
108  using namespace edm;
109  using namespace reco;
110  using namespace std;
111 
113  iEvent.getByLabel(m_clusters,cH);
114  const SiPixelClusterCollectionNew & pixelClusters = *cH.product();
115 
117  iEvent.getByLabel(m_jets,jH);
118  const edm::View<reco::Jet> & jets = *jH.product();
119 
120  CaloJetCollection selectedJets;
121  for(edm::View<reco::Jet>::const_iterator it = jets.begin() ; it != jets.end() ; it++)
122  {
123  if(it->pt() > 40 && fabs(it->eta()) < 1.6)
124  {
125  const CaloJet * ca = dynamic_cast<const CaloJet *>( &(*it));
126  if(ca ==0) abort();
127  selectedJets.push_back(*ca);
128 // std::cout << "Jet eta,phi,pt: "<< it->eta() << "," << it->phi() << "," << it->pt() << std::endl;
129  }
130  }
131 
134  iSetup.get<TkPixelCPERecord>().get(m_pixelCPE , pe );
135  pp = pe.product();
136 
138  iEvent.getByLabel(m_beamSpot,beamSpot);
139 
141  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
142  const TrackerGeometry * trackerGeometry = tracker.product();
143 
144 
145  float lengthBmodule=6.66;//cm
146  std::vector<float> zProjections;
147  for(CaloJetCollection::const_iterator jit = selectedJets.begin() ; jit != selectedJets.end() ; jit++)
148  {
149  float px=jit->px();
150  float py=jit->py();
151  float pz=jit->pz();
152  float pt=jit->pt();
153 
154  float jetZOverRho = jit->momentum().Z()/jit->momentum().Rho();
155  int minSizeY = fabs(2.*jetZOverRho)-1;
156  int maxSizeY = fabs(2.*jetZOverRho)+2;
157  if( fabs(jit->eta()) > 1.6)
158  {
159  minSizeY = 1;
160  }
161 
162  for(SiPixelClusterCollectionNew::const_iterator it = pixelClusters.begin() ; it != pixelClusters.end() ; it++) //Loop on pixel modules with clusters
163  {
164  DetId id = it->detId();
165  const edmNew::DetSet<SiPixelCluster> & detset = (*it);
166  Point3DBase<float, GlobalTag> modulepos=trackerGeometry->idToDet(id)->position();
167  float zmodule = modulepos.z() - ((modulepos.x()-beamSpot->x0())*px+(modulepos.y()-beamSpot->y0())*py)/pt * pz/pt;
168  if ((fabs(deltaPhi(jit->momentum().Phi(),modulepos.phi()))< m_maxDeltaPhi*2)&&(fabs(zmodule)<(m_maxZ+lengthBmodule/2))){
169 
170  for(size_t j = 0 ; j < detset.size() ; j ++) // Loop on pixel clusters on this module
171  {
172  const SiPixelCluster & aCluster = detset[j];
173  if(aCluster.sizeX() < m_maxSizeX && aCluster.sizeY() >= minSizeY && aCluster.sizeY() <= maxSizeY) {
174  Point3DBase<float, GlobalTag> v = trackerGeometry->idToDet(id)->surface().toGlobal(pp->localParametersV( aCluster,( *trackerGeometry->idToDetUnit(id)))[0].first) ;
175  GlobalPoint v_bs(v.x()-beamSpot->x0(),v.y()-beamSpot->y0(),v.z());
176  if(fabs(deltaPhi(jit->momentum().Phi(),v_bs.phi())) < m_maxDeltaPhi)
177  {
178  float z = v.z() - ((v.x()-beamSpot->x0())*px+(v.y()-beamSpot->y0())*py)/pt * pz/pt;
179  if(fabs(z) < m_maxZ)
180  {
181  zProjections.push_back(z);
182  }
183  }
184  } //if compatible cluster
185  } // loop on module hits
186  } // if compatible module
187  } // loop on pixel modules
188 
189  } // loop on selected jets
190  std::sort(zProjections.begin(),zProjections.end());
191 
192  std::vector<float>::iterator itCenter = zProjections.begin();
193  std::vector<float>::iterator itLeftSide = zProjections.begin();
194  std::vector<float>::iterator itRightSide = zProjections.begin();
195  std::vector<int> counts;
196  float zCluster = m_clusterLength/2.0; //cm
197  int max=0;
198  std::vector<float>::iterator left,right;
199  for(;itCenter!=zProjections.end(); itCenter++)
200  {
201 
202  while(itLeftSide != zProjections.end() && (*itCenter - *itLeftSide) > zCluster ) itLeftSide++;
203  while(itRightSide != zProjections.end() && (*itRightSide - *itCenter) < zCluster ) itRightSide++;
204 
205  int n= itRightSide-itLeftSide;
206  // std::cout << "algo :"<< *itCenter << " " << itCenter-zProjections.begin() << " dists: " << (*itCenter - *itLeftSide) << " " << (*itRightSide - *itCenter) << " count: " << n << std::endl;
207  counts.push_back(n);
208  if(n > max) {
209  max=n;
210  left=itLeftSide;
211  }
212  if(n >= max) {
213  max=n;
214  right=itRightSide;
215 // std::cout << "algo :"<< i << " " << j << " " << *itCenter << " " << itCenter-zProjections.begin() << " dists: " << (*itCenter - *itLeftSide) << " " << (*itRightSide - *itCenter) << " count: " << n << std::endl;
216  }
217  }
218 
219 
220 
221 
222  float res=0;
223  if(zProjections.size() > 0)
224  {
225  res=*(left+(right-left)/2);
226 // std::cout << "RES " << res << std::endl;
227  Vertex::Error e;
228  e(0, 0) = 0.0015 * 0.0015;
229  e(1, 1) = 0.0015 * 0.0015;
230  e(2, 2) = 1.5 * 1.5;
231  Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res);
232  Vertex thePV(p, e, 1, 1, 0);
233  std::auto_ptr<reco::VertexCollection> pOut(new reco::VertexCollection());
234  pOut->push_back(thePV);
235  iEvent.put(pOut);
236  } else
237  {
238  // std::cout << "DUMMY " << res << std::endl;
239 
241  e(0, 0) = 0.0015 * 0.0015;
242  e(1, 1) = 0.0015 * 0.0015;
243  e(2, 2) = 1.5 * 1.5;
244  Vertex::Point p(beamSpot->x(res), beamSpot->y(res), res);
245  Vertex thePV(p, e, 0, 0, 0);
246  std::auto_ptr<reco::VertexCollection> pOut(new reco::VertexCollection());
247  pOut->push_back(thePV);
248  iEvent.put(pOut);
249 
250  }
251 
252 
253 }
254 
255 
256 //define this as a plug-in
T getParameter(std::string const &) const
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Jets made from CaloTowers.
Definition: CaloJet.h:29
tuple pp
Definition: createTree.py:15
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
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
virtual void produce(edm::Event &, const edm::EventSetup &) override
float float float z
int iEvent
Definition: GenABIO.cc:243
const T & max(const T &a, const T &b)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
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
FastPrimaryVertexProducer(const edm::ParameterSet &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:55
int sizeY() const
Pixel cluster – collection of neighboring pixels above threshold.
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