CMS 3D CMS Logo

PVFitter.cc
Go to the documentation of this file.
1 
14 
20 
23 
30 
33 
34 #include "Minuit2/FCNBase.h"
35 #include "Minuit2/FunctionMinimum.h"
36 #include "Minuit2/MnMigrad.h"
37 #include "Minuit2/MnPrint.h" // Defines operator<< for cout << ierr (Dario)
38 #include "TF1.h"
39 #include "TMinuitMinimizer.h"
40 
41 #include <iostream> // Dario
42 using namespace std ; // Dario
43 // ----------------------------------------------------------------------
44 // Useful function:
45 // ----------------------------------------------------------------------
46 
47 // static char * formatTime(const std::time_t & t) {
48 // struct std::tm * ptm;
49 // ptm = gmtime(&t);
50 // static char ts[32];
51 // strftime(ts,sizeof(ts),"%Y.%m.%d %H:%M:%S %Z",ptm);
52 // return ts;
53 // }
55  edm::ConsumesCollector &&iColl)
56  : ftree_(0)
57 {
58  initialize(iConfig, iColl);
59 }
60 
63  :ftree_(0)
64 {
65  initialize(iConfig, iColl);
66 }
67 
70 {
71  //In order to make fitting ROOT histograms thread safe
72  // one must call this undocumented function
73  TMinuitMinimizer::UseStaticMinuit(false);
74  debug_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("Debug");
76  iConfig.getParameter<edm::ParameterSet>("PVFitter")
77  .getUntrackedParameter<edm::InputTag>("VertexCollection",
78  edm::InputTag("offlinePrimaryVertices")));
79  do3DFit_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("Apply3DFit");
80  //writeTxt_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("WriteAscii");
81  //outputTxt_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<std::string>("AsciiFileName");
82 
83  maxNrVertices_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<unsigned int>("maxNrStoredVertices");
84  minNrVertices_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<unsigned int>("minNrVerticesForFit");
85  minVtxNdf_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("minVertexNdf");
86  maxVtxNormChi2_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("maxVertexNormChi2");
87  minVtxTracks_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<unsigned int>("minVertexNTracks");
88  minVtxWgt_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("minVertexMeanWeight");
89  maxVtxR_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("maxVertexR");
90  maxVtxZ_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("maxVertexZ");
91  errorScale_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("errorScale");
92  sigmaCut_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("nSigmaCut");
93  fFitPerBunchCrossing=iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("FitPerBunchCrossing");
94 
95  // preset quality cut to "infinite"
96  dynamicQualityCut_ = 1.e30;
97 
98  hPVx = new TH2F("hPVx","PVx vs PVz distribution",200,-maxVtxR_, maxVtxR_, 200, -maxVtxZ_, maxVtxZ_);
99  hPVy = new TH2F("hPVy","PVy vs PVz distribution",200,-maxVtxR_, maxVtxR_, 200, -maxVtxZ_, maxVtxZ_);
100 }
101 
103 
104 }
105 
106 
108 {
109 
110 // frun = iEvent.id().run();
111 // const edm::TimeValue_t ftimestamp = iEvent.time().value();
112 // const std::time_t ftmptime = ftimestamp >> 32;
113 
114 // if (fbeginLumiOfFit == -1) freftime[0] = freftime[1] = ftmptime;
115 // if (freftime[0] == 0 || ftmptime < freftime[0]) freftime[0] = ftmptime;
116 // const char* fbeginTime = formatTime(freftime[0]);
117 // sprintf(fbeginTimeOfFit,"%s",fbeginTime);
118 
119 // if (freftime[1] == 0 || ftmptime > freftime[1]) freftime[1] = ftmptime;
120 // const char* fendTime = formatTime(freftime[1]);
121 // sprintf(fendTimeOfFit,"%s",fendTime);
122 
123 // flumi = iEvent.luminosityBlock();
124 // frunFit = frun;
125 
126 // if (fbeginLumiOfFit == -1 || fbeginLumiOfFit > flumi) fbeginLumiOfFit = flumi;
127 // if (fendLumiOfFit == -1 || fendLumiOfFit < flumi) fendLumiOfFit = flumi;
128 // std::cout << "flumi = " <<flumi<<"; fbeginLumiOfFit = " << fbeginLumiOfFit <<"; fendLumiOfFit = "<<fendLumiOfFit<<std::endl;
129 
130  //------ Primary Vertices
132  bool hasPVs = false;
133  //edm::View<reco::Vertex> vertices;
134  //const reco::VertexCollection & vertices = 0;
135 
136  if ( iEvent.getByToken(vertexToken_, PVCollection ) ) {
137  //pv = *PVCollection;
138  //vertices = *PVCollection;
139  hasPVs = true;
140  }
141  //------
142 
143  if ( hasPVs ) {
144 
145  for (reco::VertexCollection::const_iterator pv = PVCollection->begin(); pv != PVCollection->end(); ++pv ) {
146 
147 
148  //for ( size_t ipv=0; ipv != pv.size(); ++ipv ) {
149 
150  //--- vertex selection
151  if ( pv->isFake() || pv->tracksSize()==0 ) continue;
152  if ( pv->ndof() < minVtxNdf_ || (pv->ndof()+3.)/pv->tracksSize()<2*minVtxWgt_ ) continue;
153  //---
154 
155  hPVx->Fill( pv->x(), pv->z() );
156  hPVy->Fill( pv->y(), pv->z() );
157 
158  //
159  // 3D fit section
160  //
161  // apply additional quality cut
162  if ( pvQuality(*pv)>dynamicQualityCut_ ) continue;
163  // if store exceeds max. size: reduce size and apply new quality cut
164  if ( pvStore_.size()>=maxNrVertices_ ) {
165  compressStore();
166  if ( pvQuality(*pv)>dynamicQualityCut_ ) continue;
167  }
168  //
169  // copy PV to store
170  //
171  int bx = iEvent.bunchCrossing();
172  BeamSpotFitPVData pvData;
173  pvData.bunchCrossing = bx;
174  pvData.position[0] = pv->x();
175  pvData.position[1] = pv->y();
176  pvData.position[2] = pv->z();
177  pvData.posError[0] = pv->xError();
178  pvData.posError[1] = pv->yError();
179  pvData.posError[2] = pv->zError();
180  pvData.posCorr[0] = pv->covariance(0,1)/pv->xError()/pv->yError();
181  pvData.posCorr[1] = pv->covariance(0,2)/pv->xError()/pv->zError();
182  pvData.posCorr[2] = pv->covariance(1,2)/pv->yError()/pv->zError();
183  pvStore_.push_back(pvData);
184 
185  if(ftree_ != 0){
186  theBeamSpotTreeData_.run(iEvent.id().run());
190  ftree_->Fill();
191  }
192 
193  if (fFitPerBunchCrossing) bxMap_[bx].push_back(pvData);
194 
195  }
196 
197  }
198 
199 
200 
201 
202 }
203 
204 void PVFitter::setTree(TTree* tree){
205  ftree_ = tree;
207 }
208 
210 
211  using namespace ROOT::Minuit2;
212  edm::LogInfo("PVFitter") << " Number of bunch crossings: " << bxMap_.size() << std::endl;
213 
214  bool fit_ok = true;
215 
216  for ( std::map<int,std::vector<BeamSpotFitPVData> >::const_iterator pvStore = bxMap_.begin();
217  pvStore!=bxMap_.end(); ++pvStore) {
218 
219  // first set null beam spot in case
220  // fit fails
221  fbspotMap[pvStore->first] = reco::BeamSpot();
222 
223  edm::LogInfo("PVFitter") << " Number of PVs collected for PVFitter: " << (pvStore->second).size() << " in bx: " << pvStore->first << std::endl;
224 
225  if ( (pvStore->second).size() <= minNrVertices_ ) {
226  edm::LogWarning("PVFitter") << " not enough PVs, continue" << std::endl;
227  fit_ok = false;
228  continue;
229  }
230 
231  //bool fit_ok = false;
232  edm::LogInfo("PVFitter") << "Calculating beam spot with PVs ..." << std::endl;
233 
234  //
235  // LL function and fitter
236  //
237  FcnBeamSpotFitPV* fcn = new FcnBeamSpotFitPV(pvStore->second);
238  //
239  // fit parameters: positions, widths, x-y correlations, tilts in xz and yz
240  //
241  MnUserParameters upar;
242  upar.Add("x" , 0. , 0.02 , -10., 10.); // 0
243  upar.Add("y" , 0. , 0.02 , -10., 10.); // 1
244  upar.Add("z" , 0. , 0.20 , -30., 30.); // 2
245  upar.Add("ex" , 0.015, 0.01 , 0. , 10.); // 3
246  upar.Add("corrxy", 0. , 0.02 , -1. , 1. ); // 4
247  upar.Add("ey" , 0.015, 0.01 , 0. , 10.); // 5
248  upar.Add("dxdz" , 0. , 0.0002, -0.1, 0.1); // 6
249  upar.Add("dydz" , 0. , 0.0002, -0.1, 0.1); // 7
250  upar.Add("ez" , 1. , 0.1 , 0. , 30.); // 8
251  upar.Add("scale", errorScale_ , errorScale_/10.,
252  errorScale_/2., errorScale_*2.); // 9
253  MnMigrad migrad(*fcn, upar);
254 
255  //
256  // first iteration without correlations
257  //
258  upar.Fix(4);
259  upar.Fix(6);
260  upar.Fix(7);
261  upar.Fix(9);
262  FunctionMinimum ierr = migrad(0,1.);
263  if ( !ierr.IsValid() ) {
264  edm::LogInfo("PVFitter") << "3D beam spot fit failed in 1st iteration" << std::endl;
265  fit_ok = false;
266  continue;
267  }
268  //
269  // refit with harder selection on vertices
270  //
271  fcn->setLimits(upar.Value(0)-sigmaCut_*upar.Value(3),
272  upar.Value(0)+sigmaCut_*upar.Value(3),
273  upar.Value(1)-sigmaCut_*upar.Value(5),
274  upar.Value(1)+sigmaCut_*upar.Value(5),
275  upar.Value(2)-sigmaCut_*upar.Value(8),
276  upar.Value(2)+sigmaCut_*upar.Value(8));
277  ierr = migrad(0,1.);
278  if ( !ierr.IsValid() ) {
279  edm::LogInfo("PVFitter") << "3D beam spot fit failed in 2nd iteration" << std::endl;
280  fit_ok = false;
281  continue;
282  }
283  //
284  // refit with correlations
285  //
286  upar.Release(4);
287  upar.Release(6);
288  upar.Release(7);
289  ierr = migrad(0,1.);
290  if ( !ierr.IsValid() ) {
291  edm::LogInfo("PVFitter") << "3D beam spot fit failed in 3rd iteration" << std::endl;
292  fit_ok = false;
293  continue;
294  }
295  // refit with floating scale factor
296  // minuitx.ReleaseParameter(9);
297  // minuitx.Minimize();
298 
299  //minuitx.PrintResults(0,0);
300 
301  fwidthX = upar.Value(3);
302  fwidthY = upar.Value(5);
303  fwidthZ = upar.Value(8);
304  fwidthXerr = upar.Error(3);
305  fwidthYerr = upar.Error(5);
306  fwidthZerr = upar.Error(8);
307 
309  // need to get the full cov matrix
310  matrix(0,0) = pow( upar.Error(0), 2);
311  matrix(1,1) = pow( upar.Error(1), 2);
312  matrix(2,2) = pow( upar.Error(2), 2);
313  matrix(3,3) = fwidthZerr * fwidthZerr;
314  matrix(4,4) = pow( upar.Error(6), 2);
315  matrix(5,5) = pow( upar.Error(7), 2);
316  matrix(6,6) = fwidthXerr * fwidthXerr;
317 
319  upar.Value(1),
320  upar.Value(2) ),
321  fwidthZ,
322  upar.Value(6), upar.Value(7),
323  fwidthX,
324  matrix );
328 
329  fbspotMap[pvStore->first] = fbeamspot;
330  edm::LogInfo("PVFitter") << "3D PV fit done for this bunch crossing."<<std::endl;
331  //delete fcn;
332  fit_ok = fit_ok & true;
333  }
334 
335  return fit_ok;
336 }
337 
338 
340 
341  using namespace ROOT::Minuit2;
342  edm::LogInfo("PVFitter") << " Number of PVs collected for PVFitter: " << pvStore_.size() << std::endl;
343 
344  if ( pvStore_.size() <= minNrVertices_ ) return false;
345 
346  //bool fit_ok = false;
347 
348  TH1F *h1PVx = (TH1F*) hPVx->ProjectionX("h1PVx", 0, -1, "e");
349  TH1F *h1PVy = (TH1F*) hPVy->ProjectionX("h1PVy", 0, -1, "e");
350  TH1F *h1PVz = (TH1F*) hPVx->ProjectionY("h1PVz", 0, -1, "e");
351 
352  //Use our own copy for thread safety
353  TF1 gausx("localGausX","gaus");
354  TF1 gausy("localGausY","gaus");
355  TF1 gausz("localGausZ","gaus");
356 
357  h1PVx->Fit(&gausx,"QLMN0");
358  h1PVy->Fit(&gausy,"QLMN0");
359  h1PVz->Fit(&gausz,"QLMN0");
360 
361 
362  fwidthX = gausx.GetParameter(2);
363  fwidthY = gausy.GetParameter(2);
364  fwidthZ = gausz.GetParameter(2);
365  fwidthXerr = gausx.GetParError(2);
366  fwidthYerr = gausy.GetParError(2);
367  fwidthZerr = gausz.GetParError(2);
368 
369  double estX = gausx.GetParameter(1);
370  double estY = gausy.GetParameter(1);
371  double estZ = gausz.GetParameter(1);
372  double errX = fwidthX*3.;
373  double errY = fwidthY*3.;
374  double errZ = fwidthZ*3.;
375 
376  if ( ! do3DFit_ ) {
377 
379  matrix(2,2) = gausz.GetParError(1) * gausz.GetParError(1);
380  matrix(3,3) = fwidthZerr * fwidthZerr;
381  matrix(6,6) = fwidthXerr * fwidthXerr;
382 
383  fbeamspot = reco::BeamSpot( reco::BeamSpot::Point(gausx.GetParameter(1),
384  gausy.GetParameter(1),
385  gausz.GetParameter(1) ),
386  fwidthZ,
387  0., 0.,
388  fwidthX,
389  matrix );
393 
394  }
395  else { // do 3D fit
396  //
397  // LL function and fitter
398  //
400  //
401  // fit parameters: positions, widths, x-y correlations, tilts in xz and yz
402  //
403  MnUserParameters upar;
404  upar.Add("x" , estX , errX , -10. , 10. ); // 0
405  upar.Add("y" , estY , errY , -10. , 10. ); // 1
406  upar.Add("z" , estZ , errZ , -30. , 30. ); // 2
407  upar.Add("ex" , 0.015 , 0.01 , 0. , 10. ); // 3
408  upar.Add("corrxy", 0. , 0.02 , -1. , 1. ); // 4
409  upar.Add("ey" , 0.015 , 0.01 , 0. , 10. ); // 5
410  upar.Add("dxdz" , 0. , 0.0002 , -0.1 , 0.1 ); // 6
411  upar.Add("dydz" , 0. , 0.0002 , -0.1 , 0.1 ); // 7
412  upar.Add("ez" , 1. , 0.1 , 0. , 30. ); // 8
413  upar.Add("scale" , errorScale_, errorScale_/10.,errorScale_/2., errorScale_*2.); // 9
414  MnMigrad migrad(*fcn, upar);
415  //
416  // first iteration without correlations
417  //
418  migrad.Fix(4);
419  migrad.Fix(6);
420  migrad.Fix(7);
421  migrad.Fix(9);
422  FunctionMinimum ierr = migrad(0,1.);
423  if ( !ierr.IsValid() ) {
424  edm::LogWarning("PVFitter") << "3D beam spot fit failed in 1st iteration" << std::endl;
425  return false;
426  }
427  //
428  // refit with harder selection on vertices
429  //
430 
431  vector<double> results ;
432  vector<double> errors ;
433  results = ierr.UserParameters().Params() ; \
434  errors = ierr.UserParameters().Errors() ; \
435 
436  fcn->setLimits(results[0]-sigmaCut_*results[3],
437  results[0]+sigmaCut_*results[3],
438  results[1]-sigmaCut_*results[5],
439  results[1]+sigmaCut_*results[5],
440  results[2]-sigmaCut_*results[8],
441  results[2]+sigmaCut_*results[8]);
442  ierr = migrad(0,1.);
443  if ( !ierr.IsValid() ) {
444  edm::LogWarning("PVFitter") << "3D beam spot fit failed in 2nd iteration" << std::endl;
445  return false;
446  }
447  //
448  // refit with correlations
449  //
450  migrad.Release(4);
451  migrad.Release(6);
452  migrad.Release(7);
453  ierr = migrad(0,1.);
454  if ( !ierr.IsValid() ) {
455  edm::LogWarning("PVFitter") << "3D beam spot fit failed in 3rd iteration" << std::endl;
456  return false;
457  }
458  // refit with floating scale factor
459  // minuitx.ReleaseParameter(9);
460  // minuitx.Minimize();
461 
462  //minuitx.PrintResults(0,0);
463 
464  results = ierr.UserParameters().Params() ; \
465  errors = ierr.UserParameters().Errors() ; \
466 
467  fwidthX = results[3];
468  fwidthY = results[5];
469  fwidthZ = results[8];
470  fwidthXerr = errors[3];
471  fwidthYerr = errors[5];
472  fwidthZerr = errors[8];
473 
474  // check errors on widths and sigmaZ for nan
476  edm::LogWarning("PVFitter") << "3D beam spot fit returns nan in 3rd iteration" << std::endl;
477  return false;
478  }
479 
481  // need to get the full cov matrix
482  matrix(0,0) = pow( errors[0], 2);
483  matrix(1,1) = pow( errors[1], 2);
484  matrix(2,2) = pow( errors[2], 2);
485  matrix(3,3) = fwidthZerr * fwidthZerr;
486  matrix(6,6) = fwidthXerr * fwidthXerr;
487 
489  results[1],
490  results[2] ),
491  fwidthZ,
492  results[6], results[7],
493  fwidthX,
494  matrix );
498  }
499 
500  return true; //FIXME: Need to add quality test for the fit results!
501 }
502 
504 /*
505  fasciiFile << "Runnumber " << frun << std::endl;
506  fasciiFile << "BeginTimeOfFit " << fbeginTimeOfFit << std::endl;
507  fasciiFile << "EndTimeOfFit " << fendTimeOfFit << std::endl;
508  fasciiFile << "LumiRange " << fbeginLumiOfFit << " - " << fendLumiOfFit << std::endl;
509  fasciiFile << "Type " << fbeamspot.type() << std::endl;
510  fasciiFile << "X0 " << fbeamspot.x0() << std::endl;
511  fasciiFile << "Y0 " << fbeamspot.y0() << std::endl;
512  fasciiFile << "Z0 " << fbeamspot.z0() << std::endl;
513  fasciiFile << "sigmaZ0 " << fbeamspot.sigmaZ() << std::endl;
514  fasciiFile << "dxdz " << fbeamspot.dxdz() << std::endl;
515  fasciiFile << "dydz " << fbeamspot.dydz() << std::endl;
516  if (inputBeamWidth_ > 0 ) {
517  fasciiFile << "BeamWidthX " << inputBeamWidth_ << std::endl;
518  fasciiFile << "BeamWidthY " << inputBeamWidth_ << std::endl;
519  } else {
520  fasciiFile << "BeamWidthX " << fbeamspot.BeamWidthX() << std::endl;
521  fasciiFile << "BeamWidthY " << fbeamspot.BeamWidthY() << std::endl;
522  }
523 
524  for (int i = 0; i<6; ++i) {
525  fasciiFile << "Cov("<<i<<",j) ";
526  for (int j=0; j<7; ++j) {
527  fasciiFile << fbeamspot.covariance(i,j) << " ";
528  }
529  fasciiFile << std::endl;
530  }
531  // beam width error
532  if (inputBeamWidth_ > 0 ) {
533  fasciiFile << "Cov(6,j) 0 0 0 0 0 0 " << "1e-4" << std::endl;
534  } else {
535  fasciiFile << "Cov(6,j) 0 0 0 0 0 0 " << fbeamspot.covariance(6,6) << std::endl;
536  }
537  fasciiFile << "EmittanceX " << fbeamspot.emittanceX() << std::endl;
538  fasciiFile << "EmittanceY " << fbeamspot.emittanceY() << std::endl;
539  fasciiFile << "BetaStar " << fbeamspot.betaStar() << std::endl;
540 
541 */
542 }
543 
544 
545 void
547 {
548  //
549  // fill vertex qualities
550  //
551  pvQualities_.resize(pvStore_.size());
552  for ( unsigned int i=0; i<pvStore_.size(); ++i ) pvQualities_[i] = pvQuality(pvStore_[i]);
553  sort(pvQualities_.begin(),pvQualities_.end());
554  //
555  // Set new quality cut to median. This cut will be used to reduce the
556  // number of vertices in the store and also apply to all new vertices
557  // until the next reset
558  //
560  //
561  // remove all vertices failing the cut from the store
562  // (to be moved to a more efficient memory management!)
563  //
564  unsigned int iwrite(0);
565  for ( unsigned int i=0; i<pvStore_.size(); ++i ) {
566  if ( pvQuality(pvStore_[i])>dynamicQualityCut_ ) continue;
567  if ( i!=iwrite ) pvStore_[iwrite] = pvStore_[i];
568  ++iwrite;
569  }
570  pvStore_.resize(iwrite);
571  edm::LogInfo("PVFitter") << "Reduced primary vertex store size to "
572  << pvStore_.size() << " ; new dynamic quality cut = "
573  << dynamicQualityCut_ << std::endl;
574 
575 }
576 
577 double
579 {
580  //
581  // determinant of the transverse part of the PV covariance matrix
582  //
583  return
584  pv.covariance(0,0)*pv.covariance(1,1)-
585  pv.covariance(0,1)*pv.covariance(0,1);
586 }
587 
588 double
590 {
591  //
592  // determinant of the transverse part of the PV covariance matrix
593  //
594  double ex = pv.posError[0];
595  double ey = pv.posError[1];
596  return ex*ex*ey*ey*(1-pv.posCorr[0]*pv.posCorr[0]);
597 }
RunNumber_t run() const
Definition: EventID.h:39
size
Write out results.
std::vector< double > pvQualities_
Definition: PVFitter.h:205
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
T getParameter(std::string const &) const
bool debug_
Definition: PVFitter.h:145
double maxVtxNormChi2_
Definition: PVFitter.h:154
void setTree(TTree *tree)
Definition: PVFitter.cc:204
std::map< int, std::vector< BeamSpotFitPVData > > bxMap_
Definition: PVFitter.h:203
void initialize(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iColl)
Definition: PVFitter.cc:68
void run(unsigned int run)
void compressStore()
reduce size of primary vertex cache by increasing quality limit
Definition: PVFitter.cc:546
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
double minVtxWgt_
Definition: PVFitter.h:156
double minVtxNdf_
Definition: PVFitter.h:153
int bunchCrossing() const
Definition: EventBase.h:64
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:61
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
double fwidthY
Definition: PVFitter.h:183
double errorScale_
Definition: PVFitter.h:159
bool do3DFit_
Definition: PVFitter.h:146
double dynamicQualityCut_
Definition: PVFitter.h:204
double fwidthYerr
Definition: PVFitter.h:186
double fwidthZ
Definition: PVFitter.h:184
double fwidthXerr
Definition: PVFitter.h:185
void setType(BeamType type)
set beam type
Definition: BeamSpot.h:131
void bunchCrossing(unsigned int bunchCrossing)
reco::BeamSpot fbeamspot
Definition: PVFitter.h:139
void lumi(unsigned int lumi)
int iEvent
Definition: GenABIO.cc:230
std::vector< BeamSpotFitPVData > pvStore_
Definition: PVFitter.h:202
void readEvent(const edm::Event &iEvent)
Definition: PVFitter.cc:107
bool isNotFinite(T x)
Definition: isFinite.h:10
void setBeamWidthY(double v)
Definition: BeamSpot.h:109
BeamSpotTreeData theBeamSpotTreeData_
Definition: PVFitter.h:207
bool fFitPerBunchCrossing
Definition: PVFitter.h:141
TH2F * hPVy
Definition: PVFitter.h:166
def pv(vc)
Definition: MetAnalyzer.py:6
void pvData(const BeamSpotFitPVData &pvData)
TTree * ftree_
Definition: PVFitter.h:168
std::map< int, reco::BeamSpot > fbspotMap
Definition: PVFitter.h:140
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
Definition: PVFitter.h:147
double pvQuality(const reco::Vertex &pv) const
vertex quality measure
Definition: PVFitter.cc:578
double maxVtxR_
Definition: PVFitter.h:157
double maxVtxZ_
Definition: PVFitter.h:158
unsigned int maxNrVertices_
Definition: PVFitter.h:151
double fwidthZerr
Definition: PVFitter.h:187
void branch(TTree *tree)
PVFitter()
Definition: PVFitter.h:45
unsigned int minVtxTracks_
Definition: PVFitter.h:155
void fcn(int &, double *, double &, double *, int)
bool runFitter()
Definition: PVFitter.cc:339
TH2F * hPVx
Definition: PVFitter.h:166
unsigned int minNrVertices_
Definition: PVFitter.h:152
void setLimits(float xmin, float xmax, float ymin, float ymax, float zmin, float zmax)
bool runBXFitter()
Definition: PVFitter.cc:209
double sigmaCut_
Definition: PVFitter.h:160
edm::EventID id() const
Definition: EventBase.h:58
void setBeamWidthX(double v)
Definition: BeamSpot.h:108
virtual ~PVFitter()
Definition: PVFitter.cc:102
double fwidthX
Definition: PVFitter.h:182
Definition: tree.py:1
void dumpTxtFile()
Definition: PVFitter.cc:503
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40