CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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"
38 #include "TF1.h"
39 #include "TDirectory.h"
40 #include "TMinuitMinimizer.h"
41 
42 #include <iostream>
43 #include <memory>
44 
45 #include <sstream>
46 using namespace std;
47 
49  initialize(iConfig, iColl);
50 }
51 
53  initialize(iConfig, iColl);
54 }
55 
57  //In order to make fitting ROOT histograms thread safe
58  // one must call this undocumented function
59  TMinuitMinimizer::UseStaticMinuit(false);
60  debug_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("Debug");
62  iConfig.getParameter<edm::ParameterSet>("PVFitter")
63  .getUntrackedParameter<edm::InputTag>("VertexCollection", edm::InputTag("offlinePrimaryVertices")));
64  do3DFit_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("Apply3DFit");
65  //writeTxt_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("WriteAscii");
66  //outputTxt_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<std::string>("AsciiFileName");
68  iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<unsigned int>("maxNrStoredVertices");
70  iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<unsigned int>("minNrVerticesForFit");
71  minVtxNdf_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("minVertexNdf");
73  iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("maxVertexNormChi2");
75  iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<unsigned int>("minVertexNTracks");
76  minVtxWgt_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("minVertexMeanWeight");
77  maxVtxR_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("maxVertexR");
78  maxVtxZ_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("maxVertexZ");
79  errorScale_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("errorScale");
80  sigmaCut_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("nSigmaCut");
82  iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("FitPerBunchCrossing");
83  useOnlyFirstPV_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<bool>("useOnlyFirstPV");
84  minSumPt_ = iConfig.getParameter<edm::ParameterSet>("PVFitter").getUntrackedParameter<double>("minSumPt");
85 
86  // preset quality cut to "infinite"
87  dynamicQualityCut_ = 1.e30;
88 
89  hPVx = std::make_unique<TH2F>("hPVx", "PVx vs PVz distribution", 200, -maxVtxR_, maxVtxR_, 200, -maxVtxZ_, maxVtxZ_);
90  hPVy = std::make_unique<TH2F>("hPVy", "PVy vs PVz distribution", 200, -maxVtxR_, maxVtxR_, 200, -maxVtxZ_, maxVtxZ_);
91  hPVx->SetDirectory(nullptr);
92  hPVy->SetDirectory(nullptr);
93 }
94 
96 
99 
100  pvFitter.addUntracked<bool>("Debug");
101  pvFitter.addUntracked<edm::InputTag>("VertexCollection", edm::InputTag("offlinePrimaryVertices"));
102  pvFitter.addUntracked<bool>("Apply3DFit");
103  pvFitter.addUntracked<unsigned int>("maxNrStoredVertices");
104  pvFitter.addUntracked<unsigned int>("minNrVerticesForFit");
105  pvFitter.addUntracked<double>("minVertexNdf");
106  pvFitter.addUntracked<double>("maxVertexNormChi2");
107  pvFitter.addUntracked<unsigned int>("minVertexNTracks");
108  pvFitter.addUntracked<double>("minVertexMeanWeight");
109  pvFitter.addUntracked<double>("maxVertexR");
110  pvFitter.addUntracked<double>("maxVertexZ");
111  pvFitter.addUntracked<double>("errorScale");
112  pvFitter.addUntracked<double>("nSigmaCut");
113  pvFitter.addUntracked<bool>("FitPerBunchCrossing");
114  pvFitter.addUntracked<bool>("useOnlyFirstPV");
115  pvFitter.addUntracked<double>("minSumPt");
116 
117  iDesc.add<edm::ParameterSetDescription>("PVFitter", pvFitter);
118 }
119 
121  //------ Primary Vertices
123  bool hasPVs = false;
124 
125  if (iEvent.getByToken(vertexToken_, PVCollection)) {
126  hasPVs = true;
127  }
128  //------
129 
130  if (hasPVs) {
131  for (reco::VertexCollection::const_iterator pv = PVCollection->begin(); pv != PVCollection->end(); ++pv) {
132  if (useOnlyFirstPV_) {
133  if (pv != PVCollection->begin())
134  break;
135  }
136 
137  //--- vertex selection
138  if (pv->isFake() || pv->tracksSize() == 0)
139  continue;
140  if (pv->ndof() < minVtxNdf_ || (pv->ndof() + 3.) / pv->tracksSize() < 2 * minVtxWgt_)
141  continue;
142  //---
143 
144  if (pv->tracksSize() < minVtxTracks_)
145  continue;
146 
147  const auto& testTrack = pv->trackRefAt(0);
148  if (testTrack.isNull() || !testTrack.isAvailable()) {
149  edm::LogInfo("") << "Track collection not found. Skipping cut on sumPt.";
150  } else {
151  double sumPt = 0;
152  for (auto iTrack = pv->tracks_begin(); iTrack != pv->tracks_end(); ++iTrack) {
153  const auto pt = (*iTrack)->pt();
154  sumPt += pt;
155  }
156  if (sumPt < minSumPt_)
157  continue;
158  }
159 
160  hPVx->Fill(pv->x(), pv->z());
161  hPVy->Fill(pv->y(), pv->z());
162 
163  //
164  // 3D fit section
165  //
166  // apply additional quality cut
168  continue;
169  // if store exceeds max. size: reduce size and apply new quality cut
170  if (pvStore_.size() >= maxNrVertices_) {
171  compressStore();
173  continue;
174  }
175  //
176  // copy PV to store
177  //
178  int bx = iEvent.bunchCrossing();
179  BeamSpotFitPVData pvData;
180  pvData.bunchCrossing = bx;
181  pvData.position[0] = pv->x();
182  pvData.position[1] = pv->y();
183  pvData.position[2] = pv->z();
184  pvData.posError[0] = pv->xError();
185  pvData.posError[1] = pv->yError();
186  pvData.posError[2] = pv->zError();
187  pvData.posCorr[0] = pv->covariance(0, 1) / pv->xError() / pv->yError();
188  pvData.posCorr[1] = pv->covariance(0, 2) / pv->xError() / pv->zError();
189  pvData.posCorr[2] = pv->covariance(1, 2) / pv->yError() / pv->zError();
190  pvStore_.push_back(pvData);
191 
192  if (ftree_ != nullptr) {
193  theBeamSpotTreeData_.run(iEvent.id().run());
197  ftree_->Fill();
198  }
199 
201  bxMap_[bx].push_back(pvData);
202  }
203  }
204 }
205 
206 void PVFitter::setTree(TTree* tree) {
207  ftree_ = tree;
209 }
210 
212  using namespace ROOT::Minuit2;
213  edm::LogInfo("PVFitter") << " Number of bunch crossings: " << bxMap_.size() << std::endl;
214 
215  bool fit_ok = true;
216 
217  for (std::map<int, std::vector<BeamSpotFitPVData> >::const_iterator pvStore = bxMap_.begin(); pvStore != bxMap_.end();
218  ++pvStore) {
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()
224  << " in bx: " << pvStore->first << std::endl;
225 
226  if ((pvStore->second).size() <= minNrVertices_) {
227  edm::LogWarning("PVFitter") << " not enough PVs, continue" << std::endl;
228  fit_ok = false;
229  continue;
230  }
231 
232  edm::LogInfo("PVFitter") << "Calculating beam spot with PVs ..." << std::endl;
233 
234  //
235  // LL function and fitter
236  //
237  FcnBeamSpotFitPV fcn(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., errorScale_ / 2., errorScale_ * 2.); // 9
252  MnMigrad migrad(fcn, upar);
253 
254  //
255  // first iteration without correlations
256  //
257  upar.Fix(4);
258  upar.Fix(6);
259  upar.Fix(7);
260  upar.Fix(9);
261  FunctionMinimum ierr = migrad(0, 1.);
262  if (!ierr.IsValid()) {
263  edm::LogInfo("PVFitter") << "3D beam spot fit failed in 1st iteration" << std::endl;
264  fit_ok = false;
265  continue;
266  }
267  //
268  // refit with harder selection on vertices
269  //
270  fcn.setLimits(upar.Value(0) - sigmaCut_ * upar.Value(3),
271  upar.Value(0) + sigmaCut_ * upar.Value(3),
272  upar.Value(1) - sigmaCut_ * upar.Value(5),
273  upar.Value(1) + sigmaCut_ * upar.Value(5),
274  upar.Value(2) - sigmaCut_ * upar.Value(8),
275  upar.Value(2) + sigmaCut_ * upar.Value(8));
276  ierr = migrad(0, 1.);
277  if (!ierr.IsValid()) {
278  edm::LogInfo("PVFitter") << "3D beam spot fit failed in 2nd iteration" << std::endl;
279  fit_ok = false;
280  continue;
281  }
282  //
283  // refit with correlations
284  //
285  upar.Release(4);
286  upar.Release(6);
287  upar.Release(7);
288  ierr = migrad(0, 1.);
289  if (!ierr.IsValid()) {
290  edm::LogInfo("PVFitter") << "3D beam spot fit failed in 3rd iteration" << std::endl;
291  fit_ok = false;
292  continue;
293  }
294 
295  fwidthX = upar.Value(3);
296  fwidthY = upar.Value(5);
297  fwidthZ = upar.Value(8);
298  fwidthXerr = upar.Error(3);
299  fwidthYerr = upar.Error(5);
300  fwidthZerr = upar.Error(8);
301 
303  // need to get the full cov matrix
304  matrix(0, 0) = pow(upar.Error(0), 2);
305  matrix(1, 1) = pow(upar.Error(1), 2);
306  matrix(2, 2) = pow(upar.Error(2), 2);
307  matrix(3, 3) = fwidthZerr * fwidthZerr;
308  matrix(4, 4) = pow(upar.Error(6), 2);
309  matrix(5, 5) = pow(upar.Error(7), 2);
310  matrix(6, 6) = fwidthXerr * fwidthXerr;
311 
312  fbeamspot = reco::BeamSpot(reco::BeamSpot::Point(upar.Value(0), upar.Value(1), upar.Value(2)),
313  fwidthZ,
314  upar.Value(6),
315  upar.Value(7),
316  fwidthX,
317  matrix);
321 
322  fbspotMap[pvStore->first] = fbeamspot;
323  edm::LogInfo("PVFitter") << "3D PV fit done for this bunch crossing." << std::endl;
324  fit_ok = fit_ok & true;
325  }
326 
327  return fit_ok;
328 }
329 
331  using namespace ROOT::Minuit2;
332  edm::LogInfo("PVFitter") << " Number of PVs collected for PVFitter: " << pvStore_.size() << std::endl;
333 
334  if (pvStore_.size() <= minNrVertices_)
335  return false;
336 
337  //need to create a unique histogram name to avoid ROOT trying to re-use one
338  // also tell ROOT to hide its global state
339  TDirectory::TContext guard{nullptr};
340  std::ostringstream str;
341  str << this;
342  std::unique_ptr<TH1D> h1PVx{hPVx->ProjectionX((std::string("h1PVx") + str.str()).c_str(), 0, -1, "e")};
343  std::unique_ptr<TH1D> h1PVy{hPVy->ProjectionX((std::string("h1PVy") + str.str()).c_str(), 0, -1, "e")};
344  std::unique_ptr<TH1D> h1PVz{hPVx->ProjectionY((std::string("h1PVz") + str.str()).c_str(), 0, -1, "e")};
345 
346  //Use our own copy for thread safety
347  // also do not add to global list of functions
348  TF1 gausx("localGausX", "gaus", 0., 1., TF1::EAddToList::kNo);
349  TF1 gausy("localGausY", "gaus", 0., 1., TF1::EAddToList::kNo);
350  TF1 gausz("localGausZ", "gaus", 0., 1., TF1::EAddToList::kNo);
351 
352  h1PVx->Fit(&gausx, "QLMN0 SERIAL");
353  h1PVy->Fit(&gausy, "QLMN0 SERIAL");
354  h1PVz->Fit(&gausz, "QLMN0 SERIAL");
355 
356  fwidthX = gausx.GetParameter(2);
357  fwidthY = gausy.GetParameter(2);
358  fwidthZ = gausz.GetParameter(2);
359  fwidthXerr = gausx.GetParError(2);
360  fwidthYerr = gausy.GetParError(2);
361  fwidthZerr = gausz.GetParError(2);
362 
363  double estX = gausx.GetParameter(1);
364  double estY = gausy.GetParameter(1);
365  double estZ = gausz.GetParameter(1);
366  double errX = fwidthX * 3.;
367  double errY = fwidthY * 3.;
368  double errZ = fwidthZ * 3.;
369 
370  if (!do3DFit_) {
372  matrix(2, 2) = gausz.GetParError(1) * gausz.GetParError(1);
373  matrix(3, 3) = fwidthZerr * fwidthZerr;
374  matrix(6, 6) = fwidthXerr * fwidthXerr;
375 
376  fbeamspot =
377  reco::BeamSpot(reco::BeamSpot::Point(gausx.GetParameter(1), gausy.GetParameter(1), gausz.GetParameter(1)),
378  fwidthZ,
379  0.,
380  0.,
381  fwidthX,
382  matrix);
386 
387  } else { // do 3D fit
388  //
389  // LL function and fitter
390  //
392  //
393  // fit parameters: positions, widths, x-y correlations, tilts in xz and yz
394  //
395  MnUserParameters upar;
396  upar.Add("x", estX, errX, -10., 10.); // 0
397  upar.Add("y", estY, errY, -10., 10.); // 1
398  upar.Add("z", estZ, errZ, -30., 30.); // 2
399  upar.Add("ex", 0.015, 0.01, 0., 10.); // 3
400  upar.Add("corrxy", 0., 0.02, -1., 1.); // 4
401  upar.Add("ey", 0.015, 0.01, 0., 10.); // 5
402  upar.Add("dxdz", 0., 0.0002, -0.1, 0.1); // 6
403  upar.Add("dydz", 0., 0.0002, -0.1, 0.1); // 7
404  upar.Add("ez", 1., 0.1, 0., 30.); // 8
405  upar.Add("scale", errorScale_, errorScale_ / 10., errorScale_ / 2., errorScale_ * 2.); // 9
406  MnMigrad migrad(fcn, upar);
407  //
408  // first iteration without correlations
409  //
410  migrad.Fix(4);
411  migrad.Fix(6);
412  migrad.Fix(7);
413  migrad.Fix(9);
414  FunctionMinimum ierr = migrad(0, 1.);
415  if (!ierr.IsValid()) {
416  edm::LogWarning("PVFitter") << "3D beam spot fit failed in 1st iteration" << std::endl;
417  return false;
418  }
419  //
420  // refit with harder selection on vertices
421  //
422 
423  vector<double> results;
424  vector<double> errors;
425  results = ierr.UserParameters().Params();
426  errors = ierr.UserParameters().Errors();
427 
428  fcn.setLimits(results[0] - sigmaCut_ * results[3],
429  results[0] + sigmaCut_ * results[3],
430  results[1] - sigmaCut_ * results[5],
431  results[1] + sigmaCut_ * results[5],
432  results[2] - sigmaCut_ * results[8],
433  results[2] + sigmaCut_ * results[8]);
434  ierr = migrad(0, 1.);
435  if (!ierr.IsValid()) {
436  edm::LogWarning("PVFitter") << "3D beam spot fit failed in 2nd iteration" << std::endl;
437  return false;
438  }
439  //
440  // refit with correlations
441  //
442  migrad.Release(4);
443  migrad.Release(6);
444  migrad.Release(7);
445  ierr = migrad(0, 1.);
446  if (!ierr.IsValid()) {
447  edm::LogWarning("PVFitter") << "3D beam spot fit failed in 3rd iteration" << std::endl;
448  return false;
449  }
450 
451  results = ierr.UserParameters().Params();
452  errors = ierr.UserParameters().Errors();
453 
454  fwidthX = results[3];
455  fwidthY = results[5];
456  fwidthZ = results[8];
457  fwidthXerr = errors[3];
458  fwidthYerr = errors[5];
459  fwidthZerr = errors[8];
460 
461  // check errors on widths and sigmaZ for nan
463  edm::LogWarning("PVFitter") << "3D beam spot fit returns nan in 3rd iteration" << std::endl;
464  return false;
465  }
466 
468  // need to get the full cov matrix
469  matrix(0, 0) = pow(errors[0], 2);
470  matrix(1, 1) = pow(errors[1], 2);
471  matrix(2, 2) = pow(errors[2], 2);
472  matrix(3, 3) = fwidthZerr * fwidthZerr;
473  matrix(6, 6) = fwidthXerr * fwidthXerr;
474 
476  reco::BeamSpot::Point(results[0], results[1], results[2]), fwidthZ, results[6], results[7], fwidthX, matrix);
480  }
481 
482  return true;
483 }
484 
486 
488  //
489  // fill vertex qualities
490  //
491  pvQualities_.resize(pvStore_.size());
492  for (unsigned int i = 0; i < pvStore_.size(); ++i)
494  sort(pvQualities_.begin(), pvQualities_.end());
495  //
496  // Set new quality cut to median. This cut will be used to reduce the
497  // number of vertices in the store and also apply to all new vertices
498  // until the next reset
499  //
501  //
502  // remove all vertices failing the cut from the store
503  // (to be moved to a more efficient memory management!)
504  //
505  unsigned int iwrite(0);
506  for (unsigned int i = 0; i < pvStore_.size(); ++i) {
508  continue;
509  if (i != iwrite)
510  pvStore_[iwrite] = pvStore_[i];
511  ++iwrite;
512  }
513  pvStore_.resize(iwrite);
514  edm::LogInfo("PVFitter") << "Reduced primary vertex store size to " << pvStore_.size()
515  << " ; new dynamic quality cut = " << dynamicQualityCut_ << std::endl;
516 }
517 
518 double PVFitter::pvQuality(const reco::Vertex& pv) const {
519  //
520  // determinant of the transverse part of the PV covariance matrix
521  //
522  return pv.covariance(0, 0) * pv.covariance(1, 1) - pv.covariance(0, 1) * pv.covariance(0, 1);
523 }
524 
526  //
527  // determinant of the transverse part of the PV covariance matrix
528  //
529  double ex = pv.posError[0];
530  double ey = pv.posError[1];
531  return ex * ex * ey * ey * (1 - pv.posCorr[0] * pv.posCorr[0]);
532 }
RunNumber_t run() const
Definition: EventID.h:38
double minSumPt_
Definition: PVFitter.h:157
std::vector< double > pvQualities_
Definition: PVFitter.h:179
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:29
bool debug_
Definition: PVFitter.h:142
double maxVtxNormChi2_
Definition: PVFitter.h:150
void setTree(TTree *tree)
Definition: PVFitter.cc:206
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
void initialize(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iColl)
Definition: PVFitter.cc:56
void run(unsigned int run)
std::map< int, std::vector< BeamSpotFitPVData > > bxMap_
Definition: PVFitter.h:177
void compressStore()
reduce size of primary vertex cache by increasing quality limit
Definition: PVFitter.cc:487
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
double minVtxWgt_
Definition: PVFitter.h:152
bool useOnlyFirstPV_
Definition: PVFitter.h:138
dictionary results
double minVtxNdf_
Definition: PVFitter.h:149
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:27
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:148
double fwidthY
Definition: PVFitter.h:170
double errorScale_
Definition: PVFitter.h:155
bool do3DFit_
Definition: PVFitter.h:143
double dynamicQualityCut_
Definition: PVFitter.h:178
double fwidthYerr
Definition: PVFitter.h:173
double fwidthZ
Definition: PVFitter.h:171
double fwidthXerr
Definition: PVFitter.h:172
void setType(BeamType type)
set beam type
Definition: BeamSpot.h:124
void bunchCrossing(unsigned int bunchCrossing)
reco::BeamSpot fbeamspot
Definition: PVFitter.h:135
void lumi(unsigned int lumi)
int iEvent
Definition: GenABIO.cc:224
std::vector< BeamSpotFitPVData > pvStore_
Definition: PVFitter.h:176
void readEvent(const edm::Event &iEvent)
Definition: PVFitter.cc:120
void setBeamWidthY(double v)
Definition: BeamSpot.h:105
std::unique_ptr< TH2F > hPVx
Definition: PVFitter.h:161
BeamSpotTreeData theBeamSpotTreeData_
Definition: PVFitter.h:181
bool fFitPerBunchCrossing
Definition: PVFitter.h:137
void pvData(const BeamSpotFitPVData &pvData)
TTree * ftree_
Definition: PVFitter.h:164
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
Definition: PVFitter.h:144
double pvQuality(const reco::Vertex &pv) const
vertex quality measure
Definition: PVFitter.cc:518
std::map< int, reco::BeamSpot > fbspotMap
Definition: PVFitter.h:136
double maxVtxR_
Definition: PVFitter.h:153
double maxVtxZ_
Definition: PVFitter.h:154
unsigned int maxNrVertices_
Definition: PVFitter.h:147
double fwidthZerr
Definition: PVFitter.h:174
Log< level::Info, false > LogInfo
static void fillDescription(edm::ParameterSetDescription &)
Definition: PVFitter.cc:97
void branch(TTree *tree)
PVFitter()
Definition: PVFitter.h:47
unsigned int minVtxTracks_
Definition: PVFitter.h:151
std::unique_ptr< TH2F > hPVy
Definition: PVFitter.h:162
void fcn(int &, double *, double &, double *, int)
bool runFitter()
Definition: PVFitter.cc:330
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
unsigned int minNrVertices_
Definition: PVFitter.h:148
void setLimits(float xmin, float xmax, float ymin, float ymax, float zmin, float zmax)
bool runBXFitter()
Definition: PVFitter.cc:211
double sigmaCut_
Definition: PVFitter.h:156
edm::EventID id() const
Definition: EventBase.h:59
void setBeamWidthX(double v)
Definition: BeamSpot.h:104
virtual ~PVFitter()
Definition: PVFitter.cc:95
double fwidthX
Definition: PVFitter.h:169
void dumpTxtFile()
Definition: PVFitter.cc:485
Log< level::Warning, false > LogWarning
#define str(s)
tuple size
Write out results.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29