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