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 | Static Private Attributes
PrimaryVertexValidation Class Reference

#include <Alignment/OfflineValidation/plugins/PrimaryVertexValidation.cc>

Inheritance diagram for PrimaryVertexValidation:
edm::EDAnalyzer

Public Member Functions

 PrimaryVertexValidation (const edm::ParameterSet &)
 
 ~PrimaryVertexValidation ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginJob ()
 
virtual void endJob ()
 
bool hasFirstLayerPixelHits (const reco::TrackRef track)
 
bool isHit2D (const TrackingRecHit &hit) const
 
void SetVarToZero ()
 

Private Attributes

int charge_ [nMaxtracks_]
 
double chi2_ [nMaxtracks_]
 
double chi2ndof_ [nMaxtracks_]
 
float chi2normUnbiasedVertex_ [nMaxtracks_]
 
float chi2UnbiasedVertex_ [nMaxtracks_]
 
bool debug_
 
float DOFUnbiasedVertex_ [nMaxtracks_]
 
double dszFromMyVertex_ [nMaxtracks_]
 
double dxy_ [nMaxtracks_]
 
double dxyFromMyVertex_ [nMaxtracks_]
 
double dz_ [nMaxtracks_]
 
double dzFromMyVertex_ [nMaxtracks_]
 
double eta_ [nMaxtracks_]
 
std::string filename_
 
int hasRecVertex_ [nMaxtracks_]
 
int isGoodTrack_ [nMaxtracks_]
 
int Nevt_
 
int nhits1D_ [nMaxtracks_]
 
int nhits2D_ [nMaxtracks_]
 
int nhits_ [nMaxtracks_]
 
int nhitsBPIX_ [nMaxtracks_]
 
int nhitsFPIX_ [nMaxtracks_]
 
int nhitsTEC_ [nMaxtracks_]
 
int nhitsTIB_ [nMaxtracks_]
 
int nhitsTID_ [nMaxtracks_]
 
int nhitsTOB_ [nMaxtracks_]
 
int nTracks_
 
double p_ [nMaxtracks_]
 
double phi_ [nMaxtracks_]
 
double pt_ [nMaxtracks_]
 
double qoverp_ [nMaxtracks_]
 
TFile * rootFile_
 
TTree * rootTree_
 
float sumOfWeightsUnbiasedVertex_ [nMaxtracks_]
 
edm::ParameterSet theConfig
 
TrackFilterForPVFinding theTrackFilter_
 
edm::InputTag TrackCollectionTag_
 
int tracksUsedForVertexing_ [nMaxtracks_]
 
double xPCA_ [nMaxtracks_]
 
double xUnbiasedVertex_ [nMaxtracks_]
 
double yPCA_ [nMaxtracks_]
 
double yUnbiasedVertex_ [nMaxtracks_]
 
double zPCA_ [nMaxtracks_]
 
double zUnbiasedVertex_ [nMaxtracks_]
 

Static Private Attributes

static const int nMaxtracks_ = 1000
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Description: Validate alignment constants using unbiased vertex residuals

Implementation: <Notes on="" implementation>="">

Definition at line 75 of file PrimaryVertexValidation.h.

Constructor & Destructor Documentation

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

Definition at line 68 of file PrimaryVertexValidation.cc.

References debug_, filename_, edm::ParameterSet::getParameter(), and TrackCollectionTag_.

69  : theConfig(iConfig),
70  theTrackFilter_(iConfig.getParameter<edm::ParameterSet>("TkFilterParameters"))
71 {
72  //now do what ever initialization is needed
73  debug_ = iConfig.getParameter<bool> ("Debug");
74  TrackCollectionTag_ = iConfig.getParameter<edm::InputTag>("TrackCollectionTag");
75  filename_ = iConfig.getParameter<std::string>("OutputFileName");
76 }
T getParameter(std::string const &) const
TrackFilterForPVFinding theTrackFilter_
PrimaryVertexValidation::~PrimaryVertexValidation ( )

Definition at line 79 of file PrimaryVertexValidation.cc.

80 {
81  // do anything here that needs to be done at desctruction time
82  // (e.g. close files, deallocate resources etc.)
83 }

Member Function Documentation

void PrimaryVertexValidation::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
privatevirtual

Implements edm::EDAnalyzer.

Definition at line 92 of file PrimaryVertexValidation.cc.

References charge_, chi2_, chi2ndof_, chi2normUnbiasedVertex_, chi2UnbiasedVertex_, gather_cfg::cout, debug_, TransientVertex::degreesOfFreedom(), DOFUnbiasedVertex_, dszFromMyVertex_, dxy_, dxyFromMyVertex_, dz_, dzFromMyVertex_, eta_, edm::EventSetup::get(), edm::Event::getByLabel(), hasFirstLayerPixelHits(), hasRecVertex_, TransientVertex::hasTrackWeight(), i, isGoodTrack_, isHit2D(), TransientVertex::isValid(), j, kBPIX, kFPIX, LogTrace, Nevt_, nhits1D_, nhits2D_, nhits_, nhitsBPIX_, nhitsFPIX_, nhitsTEC_, nhitsTIB_, nhitsTID_, nhitsTOB_, nMaxtracks_, TransientVertex::normalisedChiSquared(), nTracks_, p_, phi_, TransientVertex::position(), pt_, qoverp_, rootTree_, SetVarToZero(), sumOfWeightsUnbiasedVertex_, StripSubdetector::TEC, theTrackFilter_, StripSubdetector::TIB, StripSubdetector::TID, StripSubdetector::TOB, TransientVertex::totalChiSquared(), ExpressReco_HICollisions_FallBack::track, TrackCollectionTag_, tracksUsedForVertexing_, TransientVertex::trackWeight(), VertexFitter< N >::vertex(), cms::Exception::what(), PV3DBase< T, PVType, FrameType >::x(), xPCA_, xUnbiasedVertex_, PV3DBase< T, PVType, FrameType >::y(), yPCA_, yUnbiasedVertex_, PV3DBase< T, PVType, FrameType >::z(), zPCA_, and zUnbiasedVertex_.

93 {
94 
95  Nevt_++;
96 
97  //=======================================================
98  // Initialize Root-tuple variables
99  //=======================================================
100 
101  SetVarToZero();
102 
103  //=======================================================
104  // Retrieve the Magnetic Field information
105  //=======================================================
106 
107  edm::ESHandle<MagneticField> theMGField;
108  iSetup.get<IdealMagneticFieldRecord>().get( theMGField );
109 
110  //=======================================================
111  // Retrieve the Tracking Geometry information
112  //=======================================================
113 
114  edm::ESHandle<GlobalTrackingGeometry> theTrackingGeometry;
115  iSetup.get<GlobalTrackingGeometryRecord>().get( theTrackingGeometry );
116 
117  //=======================================================
118  // Retrieve the Track information
119  //=======================================================
120 
121  edm::Handle<reco::TrackCollection> trackCollectionHandle;
122  iEvent.getByLabel(TrackCollectionTag_, trackCollectionHandle);
123 
124  //=======================================================
125  // Retrieve offline vartex information (only for reco)
126  //=======================================================
127 
128  /*
129  double OfflineVertexX = 0.;
130  double OfflineVertexY = 0.;
131  double OfflineVertexZ = 0.;
132  edm::Handle<reco::VertexCollection> vertices;
133  try {
134  iEvent.getByLabel("offlinePrimaryVertices", vertices);
135  } catch (...) {
136  std::cout << "No offlinePrimaryVertices found!" << std::endl;
137  }
138  if ( vertices.isValid() ) {
139  OfflineVertexX = (*vertices)[0].x();
140  OfflineVertexY = (*vertices)[0].y();
141  OfflineVertexZ = (*vertices)[0].z();
142  }
143  */
144 
145  //=======================================================
146  // Retrieve Beamspot information (only for reco)
147  //=======================================================
148 
149  /*
150  BeamSpot beamSpot;
151  edm::Handle<BeamSpot> beamSpotHandle;
152  iEvent.getByLabel("offlineBeamSpot", beamSpotHandle);
153 
154  if ( beamSpotHandle.isValid() )
155  {
156  beamSpot = *beamSpotHandle;
157 
158  } else
159  {
160  edm::LogInfo("PrimaryVertexValidation")
161  << "No beam spot available from EventSetup \n";
162  }
163 
164  double BSx0 = beamSpot.x0();
165  double BSy0 = beamSpot.y0();
166  double BSz0 = beamSpot.z0();
167 
168  if(debug_)
169  std::cout<<"Beamspot x:"<<BSx0<<" y:"<<BSy0<<" z:"<<BSz0<std::<std::endl;
170 
171  //double sigmaz = beamSpot.sigmaZ();
172  //double dxdz = beamSpot.dxdz();
173  //double BeamWidth = beamSpot.BeamWidth();
174  */
175 
176  //=======================================================
177  // Starts here ananlysis
178  //=======================================================
179 
180  if(debug_)
181  std::cout<<"PrimaryVertexValidation::analyze() looping over "<<trackCollectionHandle->size()<< "tracks." <<std::endl;
182 
183  unsigned int i = 0;
184  for(reco::TrackCollection::const_iterator track = trackCollectionHandle->begin(); track!= trackCollectionHandle->end(); ++track, ++i)
185  {
186  if ( nTracks_ >= nMaxtracks_ ) {
187  std::cout << " PrimaryVertexValidation::analyze() : Warning - Number of tracks: " << nTracks_ << " , greater than " << nMaxtracks_ << std::endl;
188  continue;
189  }
190 
191  pt_[nTracks_] = track->pt();
192  p_[nTracks_] = track->p();
193  nhits_[nTracks_] = track->numberOfValidHits();
194  eta_[nTracks_] = track->eta();
195  phi_[nTracks_] = track->phi();
196  chi2_[nTracks_] = track->chi2();
197  chi2ndof_[nTracks_] = track->normalizedChi2();
198  charge_[nTracks_] = track->charge();
199  qoverp_[nTracks_] = track->qoverp();
200  dz_[nTracks_] = track->dz();
201  dxy_[nTracks_] = track->dxy();
202  xPCA_[nTracks_] = track->vertex().x();
203  yPCA_[nTracks_] = track->vertex().y();
204  zPCA_[nTracks_] = track->vertex().z();
205 
206  //=======================================================
207  // Retrieve rechit information
208  //=======================================================
209 
210  int nRecHit1D =0;
211  int nRecHit2D =0;
212  int nhitinTIB =0;
213  int nhitinTOB =0;
214  int nhitinTID =0;
215  int nhitinTEC =0;
216  int nhitinBPIX=0;
217  int nhitinFPIX=0;
218 
219  for (trackingRecHit_iterator iHit = track->recHitsBegin(); iHit != track->recHitsEnd(); ++iHit) {
220  if((*iHit)->isValid()) {
221 
222  if (this->isHit2D(**iHit)) {++nRecHit2D;}
223  else {++nRecHit1D; }
224 
225  int type =(*iHit)->geographicalId().subdetId();
226 
227  if(type==int(StripSubdetector::TIB)){++nhitinTIB;}
228  if(type==int(StripSubdetector::TOB)){++nhitinTOB;}
229  if(type==int(StripSubdetector::TID)){++nhitinTID;}
230  if(type==int(StripSubdetector::TEC)){++nhitinTEC;}
231  if(type==int( kBPIX)){++nhitinBPIX;}
232  if(type==int( kFPIX)){++nhitinFPIX;}
233 
234  }
235  }
236 
237  nhits1D_[nTracks_] =nRecHit1D;
238  nhits2D_[nTracks_] =nRecHit2D;
239  nhitsBPIX_[nTracks_] =nhitinBPIX;
240  nhitsFPIX_[nTracks_] =nhitinFPIX;
241  nhitsTIB_[nTracks_] =nhitinTIB;
242  nhitsTID_[nTracks_] =nhitinTID;
243  nhitsTOB_[nTracks_] =nhitinTOB;
244  nhitsTEC_[nTracks_] =nhitinTEC;
245 
246  //=======================================================
247  // Good tracks for vertexing selection
248  //=======================================================
249 
250  reco::TrackRef trackref(trackCollectionHandle,i);
251  bool hasTheProbeFirstPixelLayerHit = false;
252  hasTheProbeFirstPixelLayerHit = this->hasFirstLayerPixelHits(trackref);
253  reco::TransientTrack theTTRef = reco::TransientTrack(trackref, &*theMGField, theTrackingGeometry );
254  if (theTrackFilter_(theTTRef)&&hasTheProbeFirstPixelLayerHit){
256  }
257 
258  //=======================================================
259  // Fit unbiased vertex
260  //=======================================================
261 
262  std::vector<reco::TransientTrack> transientTracks;
263  for(size_t j = 0; j < trackCollectionHandle->size(); j++)
264  {
265  reco::TrackRef tk(trackCollectionHandle, j);
266  if( tk == trackref ) continue;
267  bool hasTheTagFirstPixelLayerHit = false;
268  hasTheTagFirstPixelLayerHit = this->hasFirstLayerPixelHits(tk);
269  reco::TransientTrack theTT = reco::TransientTrack(tk, &*theMGField, theTrackingGeometry );
270  if (theTrackFilter_(theTT)&&hasTheTagFirstPixelLayerHit){
271  transientTracks.push_back(theTT);
272  }
273  }
274 
275  if(transientTracks.size() > 2){
276 
277  if(debug_)
278  std::cout <<"PrimaryVertexValidation::analyze() :Transient Track Collection size: "<<transientTracks.size()<<std::endl;
279 
280  try{
281 
283  TransientVertex theFittedVertex = fitter->vertex(transientTracks);
284 
285  if(theFittedVertex.isValid ()){
286 
287  if(theFittedVertex.hasTrackWeight()){
288  for(size_t rtracks= 0; rtracks < transientTracks.size(); rtracks++){
289  sumOfWeightsUnbiasedVertex_[nTracks_] += theFittedVertex.trackWeight(transientTracks[rtracks]);
290  }
291  }
292 
293  const math::XYZPoint myVertex(theFittedVertex.position().x(),theFittedVertex.position().y(),theFittedVertex.position().z());
294  hasRecVertex_[nTracks_] = 1;
295  xUnbiasedVertex_[nTracks_] = theFittedVertex.position().x();
296  yUnbiasedVertex_[nTracks_] = theFittedVertex.position().y();
297  zUnbiasedVertex_[nTracks_] = theFittedVertex.position().z();
299  chi2UnbiasedVertex_[nTracks_] = theFittedVertex.totalChiSquared();
300  DOFUnbiasedVertex_[nTracks_] = theFittedVertex.degreesOfFreedom();
301  tracksUsedForVertexing_[nTracks_] = transientTracks.size();
302  dxyFromMyVertex_[nTracks_] = track->dxy(myVertex);
303  dzFromMyVertex_[nTracks_] = track->dz(myVertex);
304  dszFromMyVertex_[nTracks_] = track->dsz(myVertex);
305 
306  if(debug_){
307  std::cout<<"PrimaryVertexValidation::analyze() :myVertex.x()= "<<myVertex.x()<<" myVertex.y()= "<<myVertex.y()<<" theFittedVertex.z()= "<<myVertex.z()<<std::endl;
308  std::cout<<"PrimaryVertexValidation::analyze() : track->dz(myVertex)= "<<track->dz(myVertex)<<std::endl;
309  std::cout<<"PrimaryVertexValidation::analyze() : zPCA -myVertex.z() = "<<(track->vertex().z() -myVertex.z() )<<std::endl;
310  }// ends if debug_
311  } // ends if the fitted vertex is Valid
312 
313  } catch ( cms::Exception& er ) {
314  LogTrace("PrimaryVertexValidation::analyze RECO")<<"caught std::exception "<<er.what()<<std::endl;
315  }
316  } //ends if transientTracks.size() > 2
317 
318  else {
319  std::cout << "PrimaryVertexValidation::analyze() :Not enough tracks to make a vertex. Returns no vertex info" << std::endl;
320  }
321 
322  ++nTracks_;
323 
324  if(debug_)
325  std::cout<< "Track "<<i<<" : pT = "<<track->pt()<<std::endl;
326 
327  }// for loop on tracks
328 
329  rootTree_->Fill();
330 }
virtual char const * what() const
Definition: Exception.cc:97
type
Definition: HCALResponse.h:22
int i
Definition: DBlmapReader.cc:9
float totalChiSquared() const
T y() const
Definition: PV3DBase.h:57
float DOFUnbiasedVertex_[nMaxtracks_]
double dszFromMyVertex_[nMaxtracks_]
double dxyFromMyVertex_[nMaxtracks_]
bool hasFirstLayerPixelHits(const reco::TrackRef track)
TrackFilterForPVFinding theTrackFilter_
const int kFPIX
float normalisedChiSquared() const
int tracksUsedForVertexing_[nMaxtracks_]
virtual CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const =0
float degreesOfFreedom() const
double yUnbiasedVertex_[nMaxtracks_]
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:58
int j
Definition: DBlmapReader.cc:9
double zUnbiasedVertex_[nMaxtracks_]
float chi2normUnbiasedVertex_[nMaxtracks_]
float trackWeight(const reco::TransientTrack &track) const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
#define LogTrace(id)
bool hasTrackWeight() const
bool isHit2D(const TrackingRecHit &hit) const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
const T & get() const
Definition: EventSetup.h:55
double xUnbiasedVertex_[nMaxtracks_]
double dzFromMyVertex_[nMaxtracks_]
tuple cout
Definition: gather_cfg.py:41
float sumOfWeightsUnbiasedVertex_[nMaxtracks_]
T x() const
Definition: PV3DBase.h:56
bool isValid() const
float chi2UnbiasedVertex_[nMaxtracks_]
const int kBPIX
void PrimaryVertexValidation::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 382 of file PrimaryVertexValidation.cc.

References charge_, chi2_, chi2ndof_, chi2normUnbiasedVertex_, chi2UnbiasedVertex_, DOFUnbiasedVertex_, dszFromMyVertex_, dxy_, dxyFromMyVertex_, dz_, dzFromMyVertex_, eta_, filename_, hasRecVertex_, isGoodTrack_, Nevt_, nhits1D_, nhits2D_, nhits_, nhitsBPIX_, nhitsFPIX_, nhitsTEC_, nhitsTIB_, nhitsTID_, nhitsTOB_, nTracks_, p_, phi_, pt_, qoverp_, rootFile_, rootTree_, sumOfWeightsUnbiasedVertex_, tracksUsedForVertexing_, xPCA_, xUnbiasedVertex_, yPCA_, yUnbiasedVertex_, zPCA_, and zUnbiasedVertex_.

383 {
384  edm::LogInfo("beginJob") << "Begin Job" << std::endl;
385  // Define TTree for output
386 
387  Nevt_ = 0;
388 
389  rootFile_ = new TFile(filename_.c_str(),"recreate");
390  rootTree_ = new TTree("tree","PV Validation tree");
391 
392  // Track Paramters
393  rootTree_->Branch("nTracks",&nTracks_,"nTracks/I");
394  rootTree_->Branch("pt",&pt_,"pt[nTracks]/D");
395  rootTree_->Branch("p",&p_,"p[nTracks]/D");
396  rootTree_->Branch("nhits",&nhits_,"nhits[nTracks]/I");
397  rootTree_->Branch("nhits1D",&nhits1D_,"nhits1D[nTracks]/I");
398  rootTree_->Branch("nhits2D",&nhits2D_,"nhits2D[nTracks]/I");
399  rootTree_->Branch("nhitsBPIX",&nhitsBPIX_,"nhitsBPIX[nTracks]/I");
400  rootTree_->Branch("nhitsFPIX",&nhitsFPIX_,"nhitsFPIX[nTracks]/I");
401  rootTree_->Branch("nhitsTIB",&nhitsTIB_,"nhitsTIB[nTracks]/I");
402  rootTree_->Branch("nhitsTID",&nhitsTID_,"nhitsTID[nTracks]/I");
403  rootTree_->Branch("nhitsTOB",&nhitsTOB_,"nhitsTOB[nTracks]/I");
404  rootTree_->Branch("nhitsTEC",&nhitsTEC_,"nhitsTEC[nTracks]/I");
405  rootTree_->Branch("eta",&eta_,"eta[nTracks]/D");
406  rootTree_->Branch("phi",&phi_,"phi[nTracks]/D");
407  rootTree_->Branch("chi2",&chi2_,"chi2[nTracks]/D");
408  rootTree_->Branch("chi2ndof",&chi2ndof_,"chi2ndof[nTracks]/D");
409  rootTree_->Branch("charge",&charge_,"charge[nTracks]/I");
410  rootTree_->Branch("qoverp",&qoverp_,"qoverp[nTracks]/D");
411  rootTree_->Branch("dz",&dz_,"dz[nTracks]/D");
412  rootTree_->Branch("dxy",&dxy_,"dxy[nTracks]/D");
413  rootTree_->Branch("xPCA",&xPCA_,"xPCA[nTracks]/D");
414  rootTree_->Branch("yPCA",&yPCA_,"yPCA[nTracks]/D");
415  rootTree_->Branch("zPCA",&zPCA_,"zPCA[nTracks]/D");
416  rootTree_->Branch("xUnbiasedVertex",&xUnbiasedVertex_,"xUnbiasedVertex[nTracks]/D");
417  rootTree_->Branch("yUnbiasedVertex",&yUnbiasedVertex_,"yUnbiasedVertex[nTracks]/D");
418  rootTree_->Branch("zUnbiasedVertex",&zUnbiasedVertex_,"zUnbiasedVertex[nTracks]/D");
419  rootTree_->Branch("chi2normUnbiasedVertex",&chi2normUnbiasedVertex_,"chi2normUnbiasedVertex[nTracks]/F");
420  rootTree_->Branch("chi2UnbiasedVertex",&chi2UnbiasedVertex_,"chi2UnbiasedVertex[nTracks]/F");
421  rootTree_->Branch("DOFUnbiasedVertex",&DOFUnbiasedVertex_," DOFUnbiasedVertex[nTracks]/F");
422  rootTree_->Branch("sumOfWeightsUnbiasedVertex",&sumOfWeightsUnbiasedVertex_,"sumOfWeightsUnbiasedVertex[nTracks]/F");
423  rootTree_->Branch("tracksUsedForVertexing",&tracksUsedForVertexing_,"tracksUsedForVertexing[nTracks]/I");
424  rootTree_->Branch("dxyFromMyVertex",&dxyFromMyVertex_,"dxyFromMyVertex[nTracks]/D");
425  rootTree_->Branch("dzFromMyVertex",&dzFromMyVertex_,"dzFromMyVertex[nTracks]/D");
426  rootTree_->Branch("dszFromMyVertex",&dszFromMyVertex_,"dszFromMyVertex[nTracks]/D");
427  rootTree_->Branch("hasRecVertex",&hasRecVertex_,"hasRecVertex[nTracks]/I");
428  rootTree_->Branch("isGoodTrack",&isGoodTrack_,"isGoodTrack[nTracks]/I");
429 }
float DOFUnbiasedVertex_[nMaxtracks_]
double dszFromMyVertex_[nMaxtracks_]
double dxyFromMyVertex_[nMaxtracks_]
int tracksUsedForVertexing_[nMaxtracks_]
double yUnbiasedVertex_[nMaxtracks_]
double zUnbiasedVertex_[nMaxtracks_]
float chi2normUnbiasedVertex_[nMaxtracks_]
double xUnbiasedVertex_[nMaxtracks_]
double dzFromMyVertex_[nMaxtracks_]
float sumOfWeightsUnbiasedVertex_[nMaxtracks_]
float chi2UnbiasedVertex_[nMaxtracks_]
void PrimaryVertexValidation::endJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 432 of file PrimaryVertexValidation.cc.

References gather_cfg::cout, Nevt_, and rootFile_.

433 {
434 
435  std::cout<<"######################################"<<std::endl;
436  std::cout<<"Number of analyzed events: "<<Nevt_<<std::endl;
437  std::cout<<"######################################"<<std::endl;
438 
439  if ( rootFile_ ) {
440  rootFile_->Write();
441  rootFile_->Close();
442  }
443 }
tuple cout
Definition: gather_cfg.py:41
bool PrimaryVertexValidation::hasFirstLayerPixelHits ( const reco::TrackRef  track)
private

Definition at line 363 of file PrimaryVertexValidation.cc.

References reco::HitPattern::getHitPattern(), reco::HitPattern::getLayer(), i, reco::HitPattern::numberOfHits(), L1TEmulatorMonitor_cff::p, listBenchmarks::pattern, reco::HitPattern::pixelBarrelHitFilter(), reco::HitPattern::pixelEndcapHitFilter(), and reco::HitPattern::validHitFilter().

Referenced by analyze().

364 {
365  bool accepted = false;
366  // hit pattern of the track
367  const reco::HitPattern& p = track->hitPattern();
368  for (int i=0; i<p.numberOfHits(); i++) {
369  uint32_t pattern = p.getHitPattern(i);
370  if (p.pixelBarrelHitFilter(pattern) || p.pixelEndcapHitFilter(pattern) ) {
371  if (p.getLayer(pattern) == 1) {
372  if (p.validHitFilter(pattern)) {
373  accepted = true;
374  }
375  }
376  }
377  }
378  return accepted;
379 }
int i
Definition: DBlmapReader.cc:9
uint32_t getLayer(uint32_t pattern) const
Definition: HitPattern.cc:219
bool pixelEndcapHitFilter(uint32_t pattern) const
Definition: HitPattern.cc:153
int numberOfHits() const
Definition: HitPattern.cc:312
bool pixelBarrelHitFilter(uint32_t pattern) const
Definition: HitPattern.cc:146
bool validHitFilter(uint32_t pattern) const
Definition: HitPattern.cc:264
uint32_t getHitPattern(int position) const
Definition: HitPattern.cc:86
bool PrimaryVertexValidation::isHit2D ( const TrackingRecHit hit) const
private

Definition at line 333 of file PrimaryVertexValidation.cc.

References TrackingRecHit::dimension(), TrackingRecHit::geographicalId(), kBPIX, kFPIX, and DetId::Tracker.

Referenced by analyze().

334 {
335  if (hit.dimension() < 2) {
336  return false; // some (muon...) stuff really has RecHit1D
337  } else {
338  const DetId detId(hit.geographicalId());
339  if (detId.det() == DetId::Tracker) {
340  if (detId.subdetId() == kBPIX || detId.subdetId() == kFPIX) {
341  return true; // pixel is always 2D
342  } else { // should be SiStrip now
343  if (dynamic_cast<const SiStripRecHit2D*>(&hit)) return false; // normal hit
344  else if (dynamic_cast<const SiStripMatchedRecHit2D*>(&hit)) return true; // matched is 2D
345  else if (dynamic_cast<const ProjectedSiStripRecHit2D*>(&hit)) return false; // crazy hit...
346  else {
347  edm::LogError("UnkownType") << "@SUB=AlignmentTrackSelector::isHit2D"
348  << "Tracker hit not in pixel and neither SiStripRecHit2D nor "
349  << "SiStripMatchedRecHit2D nor ProjectedSiStripRecHit2D.";
350  return false;
351  }
352  }
353  } else { // not tracker??
354  edm::LogWarning("DetectorMismatch") << "@SUB=AlignmentTrackSelector::isHit2D"
355  << "Hit not in tracker with 'official' dimension >=2.";
356  return true; // dimension() >= 2 so accept that...
357  }
358  }
359  // never reached...
360 }
virtual int dimension() const =0
const int kFPIX
Definition: DetId.h:20
DetId geographicalId() const
const int kBPIX
void PrimaryVertexValidation::SetVarToZero ( )
private

Definition at line 445 of file PrimaryVertexValidation.cc.

References charge_, chi2_, chi2ndof_, chi2normUnbiasedVertex_, chi2UnbiasedVertex_, DOFUnbiasedVertex_, dszFromMyVertex_, dxy_, dxyFromMyVertex_, dz_, dzFromMyVertex_, eta_, hasRecVertex_, i, isGoodTrack_, nhits1D_, nhits2D_, nhits_, nhitsBPIX_, nhitsFPIX_, nhitsTEC_, nhitsTIB_, nhitsTID_, nhitsTOB_, nMaxtracks_, nTracks_, p_, phi_, pt_, qoverp_, sumOfWeightsUnbiasedVertex_, tracksUsedForVertexing_, xPCA_, xUnbiasedVertex_, yPCA_, yUnbiasedVertex_, zPCA_, and zUnbiasedVertex_.

Referenced by analyze().

445  {
446 
447  nTracks_ = 0;
448  for ( int i=0; i<nMaxtracks_; ++i ) {
449  pt_[i] = 0;
450  p_[i] = 0;
451  nhits_[i] = 0;
452  nhits1D_[i] = 0;
453  nhits2D_[i] = 0;
454  nhitsBPIX_[i] = 0;
455  nhitsFPIX_[i] = 0;
456  nhitsTIB_[i] = 0;
457  nhitsTID_[i] = 0;
458  nhitsTOB_[i] = 0;
459  nhitsTEC_[i] = 0;
460  eta_[i] = 0;
461  phi_[i] = 0;
462  chi2_[i] = 0;
463  chi2ndof_[i] = 0;
464  charge_[i] = 0;
465  qoverp_[i] = 0;
466  dz_[i] = 0;
467  dxy_[i] = 0;
468  xPCA_[i] = 0;
469  yPCA_[i] = 0;
470  zPCA_[i] = 0;
471  xUnbiasedVertex_[i] =0;
472  yUnbiasedVertex_[i] =0;
473  zUnbiasedVertex_[i] =0;
479  dxyFromMyVertex_[i]=0;
480  dzFromMyVertex_[i]=0;
481  dszFromMyVertex_[i]=0;
482  hasRecVertex_[i] = 0;
483  isGoodTrack_[i] = 0;
484  }
485 }
int i
Definition: DBlmapReader.cc:9
float DOFUnbiasedVertex_[nMaxtracks_]
double dszFromMyVertex_[nMaxtracks_]
double dxyFromMyVertex_[nMaxtracks_]
int tracksUsedForVertexing_[nMaxtracks_]
double yUnbiasedVertex_[nMaxtracks_]
double zUnbiasedVertex_[nMaxtracks_]
float chi2normUnbiasedVertex_[nMaxtracks_]
double xUnbiasedVertex_[nMaxtracks_]
double dzFromMyVertex_[nMaxtracks_]
float sumOfWeightsUnbiasedVertex_[nMaxtracks_]
float chi2UnbiasedVertex_[nMaxtracks_]

Member Data Documentation

int PrimaryVertexValidation::charge_[nMaxtracks_]
private

Definition at line 121 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::chi2_[nMaxtracks_]
private

Definition at line 119 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::chi2ndof_[nMaxtracks_]
private

Definition at line 120 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

float PrimaryVertexValidation::chi2normUnbiasedVertex_[nMaxtracks_]
private

Definition at line 131 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

float PrimaryVertexValidation::chi2UnbiasedVertex_[nMaxtracks_]
private

Definition at line 132 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

bool PrimaryVertexValidation::debug_
private

Definition at line 91 of file PrimaryVertexValidation.h.

Referenced by analyze(), and PrimaryVertexValidation().

float PrimaryVertexValidation::DOFUnbiasedVertex_[nMaxtracks_]
private

Definition at line 133 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::dszFromMyVertex_[nMaxtracks_]
private

Definition at line 138 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::dxy_[nMaxtracks_]
private

Definition at line 124 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::dxyFromMyVertex_[nMaxtracks_]
private

Definition at line 136 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::dz_[nMaxtracks_]
private

Definition at line 123 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::dzFromMyVertex_[nMaxtracks_]
private

Definition at line 137 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::eta_[nMaxtracks_]
private

Definition at line 117 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

std::string PrimaryVertexValidation::filename_
private

Definition at line 96 of file PrimaryVertexValidation.h.

Referenced by beginJob(), and PrimaryVertexValidation().

int PrimaryVertexValidation::hasRecVertex_[nMaxtracks_]
private

Definition at line 139 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::isGoodTrack_[nMaxtracks_]
private

Definition at line 140 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::Nevt_
private

Definition at line 92 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and endJob().

int PrimaryVertexValidation::nhits1D_[nMaxtracks_]
private

Definition at line 109 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::nhits2D_[nMaxtracks_]
private

Definition at line 110 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::nhits_[nMaxtracks_]
private

Definition at line 108 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::nhitsBPIX_[nMaxtracks_]
private

Definition at line 111 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::nhitsFPIX_[nMaxtracks_]
private

Definition at line 112 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::nhitsTEC_[nMaxtracks_]
private

Definition at line 116 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::nhitsTIB_[nMaxtracks_]
private

Definition at line 113 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::nhitsTID_[nMaxtracks_]
private

Definition at line 114 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

int PrimaryVertexValidation::nhitsTOB_[nMaxtracks_]
private

Definition at line 115 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

const int PrimaryVertexValidation::nMaxtracks_ = 1000
staticprivate

Definition at line 104 of file PrimaryVertexValidation.h.

Referenced by analyze(), and SetVarToZero().

int PrimaryVertexValidation::nTracks_
private

Definition at line 105 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::p_[nMaxtracks_]
private

Definition at line 107 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::phi_[nMaxtracks_]
private

Definition at line 118 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::pt_[nMaxtracks_]
private

Definition at line 106 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::qoverp_[nMaxtracks_]
private

Definition at line 122 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

TFile* PrimaryVertexValidation::rootFile_
private

Definition at line 97 of file PrimaryVertexValidation.h.

Referenced by beginJob(), and endJob().

TTree* PrimaryVertexValidation::rootTree_
private

Definition at line 98 of file PrimaryVertexValidation.h.

Referenced by analyze(), and beginJob().

float PrimaryVertexValidation::sumOfWeightsUnbiasedVertex_[nMaxtracks_]
private

Definition at line 134 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

edm::ParameterSet PrimaryVertexValidation::theConfig
private

Definition at line 89 of file PrimaryVertexValidation.h.

TrackFilterForPVFinding PrimaryVertexValidation::theTrackFilter_
private

Definition at line 93 of file PrimaryVertexValidation.h.

Referenced by analyze().

edm::InputTag PrimaryVertexValidation::TrackCollectionTag_
private

Definition at line 90 of file PrimaryVertexValidation.h.

Referenced by analyze(), and PrimaryVertexValidation().

int PrimaryVertexValidation::tracksUsedForVertexing_[nMaxtracks_]
private

Definition at line 135 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::xPCA_[nMaxtracks_]
private

Definition at line 125 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::xUnbiasedVertex_[nMaxtracks_]
private

Definition at line 128 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::yPCA_[nMaxtracks_]
private

Definition at line 126 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::yUnbiasedVertex_[nMaxtracks_]
private

Definition at line 129 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::zPCA_[nMaxtracks_]
private

Definition at line 127 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().

double PrimaryVertexValidation::zUnbiasedVertex_[nMaxtracks_]
private

Definition at line 130 of file PrimaryVertexValidation.h.

Referenced by analyze(), beginJob(), and SetVarToZero().