CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
HDQMfitUtilities Class Reference

#include <HDQMfitUtilities.h>

Public Member Functions

double doGaussFit (MonitorElement *ME)
 
double doGaussFit (TH1F *)
 
double doLanGaussFit (MonitorElement *ME)
 
double doLanGaussFit (TH1F *)
 
double getFitChi ()
 
int getFitnDof ()
 
double getGaussPar (std::string s)
 
double getGaussParErr (std::string s)
 
double getLanGaussConv (std::string s)
 
double getLanGaussPar (std::string s)
 
double getLanGaussParErr (std::string s)
 
 HDQMfitUtilities ()
 
void init ()
 
 ~HDQMfitUtilities ()
 

Private Attributes

double chi2GausS
 
double epGausS [3]
 
double epLanGausS [4]
 
TF1 * gausFit
 
TF1 * langausFit
 
int32_t nDofGausS
 
double pGausS [3]
 
double pLanConv [2]
 
double pLanGausS [4]
 

Detailed Description

Definition at line 28 of file HDQMfitUtilities.h.

Constructor & Destructor Documentation

HDQMfitUtilities::HDQMfitUtilities ( )

Definition at line 173 of file HDQMfitUtilities.cc.

References init().

173 : langausFit(nullptr), gausFit(nullptr) { init(); }
HDQMfitUtilities::~HDQMfitUtilities ( )

Definition at line 202 of file HDQMfitUtilities.cc.

References gausFit, and langausFit.

202  {
203  if (langausFit != nullptr)
204  delete langausFit;
205  if (gausFit != nullptr)
206  delete gausFit;
207 }

Member Function Documentation

double HDQMfitUtilities::doGaussFit ( MonitorElement ME)
inline
double HDQMfitUtilities::doGaussFit ( TH1F *  htoFit)

Definition at line 313 of file HDQMfitUtilities.cc.

References chi2GausS, epGausS, gausFit, HDQMUtil::Gauss(), mps_fire::i, init(), langausFit, nDofGausS, pGausS, pfDeepBoostedJetPreprocessParams_cfi::sv, TrackerOfflineValidation_Dqm_cff::xmax, and Phase2TrackerMonitorDigi_cff::ymax.

313  {
314  init();
315  // if (htoFit->GetEntries()!=0) {
316  if (htoFit->Integral() != 0) {
317  // Setting fit range and start values
318  double fr[2];
319  double sv[3], pllo[3], plhi[3];
320  fr[0] = htoFit->GetMean() - 5 * htoFit->GetRMS();
321  fr[1] = htoFit->GetMean() + 5 * htoFit->GetRMS();
322 
323  int32_t imax = htoFit->GetMaximumBin();
324  double xmax = htoFit->GetBinCenter(imax);
325  double ymax = htoFit->GetBinContent(imax);
326  int32_t i[2];
327  int32_t iArea[2];
328 
329  i[0] = htoFit->GetXaxis()->FindBin(fr[0]);
330  i[1] = htoFit->GetXaxis()->FindBin(fr[1]);
331 
332  iArea[0] = htoFit->GetXaxis()->FindBin(fr[0]);
333  iArea[1] = htoFit->GetXaxis()->FindBin(fr[1]);
334  double AreaFWHM = htoFit->Integral(iArea[0], iArea[1], "width");
335 
336  sv[2] = AreaFWHM / (4 * ymax);
337  sv[1] = xmax;
338  sv[0] = htoFit->Integral(i[0], i[1], "width");
339 
340  plhi[0] = 1000000.0;
341  plhi[1] = 10.0;
342  plhi[2] = 10.;
343  pllo[0] = 1.5;
344  pllo[1] = 0.1;
345  pllo[2] = 0.3;
346  Char_t FunName[100];
347  sprintf(FunName, "FitfcnLG_%s", htoFit->GetName());
348  gausFit = new TF1(FunName, HDQMUtil::Gauss, fr[0], fr[1], 3);
349  gausFit->SetParameters(sv);
350  gausFit->SetParNames("Constant", "GaussPeak", "Sigma");
351 
352  for (int32_t i = 0; i < 3; i++) {
353  gausFit->SetParLimits(i, pllo[i], plhi[i]);
354  }
355 
356  try {
357  htoFit->Fit(gausFit, "R0");
358 
359  gausFit->SetRange(fr[0], fr[1]);
360  gausFit->GetParameters(pGausS);
361  std::memcpy((void *)epGausS, (void *)gausFit->GetParErrors(), 3 * sizeof(double));
362 
363  chi2GausS = langausFit->GetChisquare(); // obtain chi^2
364  nDofGausS = langausFit->GetNDF(); // obtain ndf
365  } catch (cms::Exception &iException) {
366  edm::LogError("fitUtility") << "problem in fitting " << htoFit->GetTitle()
367  << " \n\tDefault values of the parameters will be used";
368  pGausS[0] = -9999;
369  pGausS[1] = -9999;
370  pGausS[2] = -9999;
371  epGausS[0] = -9999;
372  epGausS[1] = -9999;
373  epGausS[2] = -9999;
374  chi2GausS = -9999;
375  nDofGausS = -9999;
376  }
377 
378  } else {
379  pGausS[0] = -9999;
380  pGausS[1] = -9999;
381  pGausS[2] = -9999;
382  epGausS[0] = -9999;
383  epGausS[1] = -9999;
384  epGausS[2] = -9999;
385  chi2GausS = -9999;
386  nDofGausS = -9999;
387  }
388 
389  return htoFit->GetEntries();
390 }
double Gauss(double *x, double *par)
double HDQMfitUtilities::doLanGaussFit ( MonitorElement ME)
inline

Definition at line 34 of file HDQMfitUtilities.h.

References doLanGaussFit(), and MonitorElement::getTH1F().

Referenced by doLanGaussFit(), and SiStripDQMHistoryHelper::setDBValuesForLandau().

34 { return doLanGaussFit(ME->getTH1F()); }
TH1F * getTH1F() const
double doLanGaussFit(MonitorElement *ME)
double HDQMfitUtilities::doLanGaussFit ( TH1F *  htoFit)

Definition at line 210 of file HDQMfitUtilities.cc.

References chi2GausS, epLanGausS, mps_fire::i, init(), HDQMUtil::langaufun(), HDQMUtil::langaupro(), langausFit, nDofGausS, pLanConv, pLanGausS, pfDeepBoostedJetPreprocessParams_cfi::sv, TrackerOfflineValidation_Dqm_cff::xmax, and Phase2TrackerMonitorDigi_cff::ymax.

210  {
211  init();
212 
213  // if (htoFit->GetEntries()!=0) {
214  // Check for the entries excluding over/underflows
215  if (htoFit->Integral() != 0) {
216  edm::LogInfo("fitUtility") << "Fitting " << htoFit->GetTitle() << std::endl;
217  // Setting fit range and start values
218  double fr[2];
219  double sv[4], pllo[4], plhi[4];
220  fr[0] = 0.5 * htoFit->GetMean();
221  fr[1] = 3.0 * htoFit->GetMean();
222 
223  // (EM) parameters setting good for signal only
224  int32_t imax = htoFit->GetMaximumBin();
225  double xmax = htoFit->GetBinCenter(imax);
226  double ymax = htoFit->GetBinContent(imax);
227  int32_t i[2];
228  int32_t iArea[2];
229 
230  i[0] = htoFit->GetXaxis()->FindBin(fr[0]);
231  i[1] = htoFit->GetXaxis()->FindBin(fr[1]);
232 
233  iArea[0] = htoFit->GetXaxis()->FindBin(fr[0]);
234  iArea[1] = htoFit->GetXaxis()->FindBin(fr[1]);
235  double AreaFWHM = htoFit->Integral(iArea[0], iArea[1], "width");
236 
237  sv[1] = xmax;
238  sv[2] = htoFit->Integral(i[0], i[1], "width");
239  sv[3] = AreaFWHM / (4 * ymax);
240  sv[0] = sv[3];
241 
242  plhi[0] = 25.0;
243  plhi[1] = 200.0;
244  plhi[2] = 1000000.0;
245  plhi[3] = 50.0;
246  pllo[0] = 1.5;
247  pllo[1] = 10.0;
248  pllo[2] = 1.0;
249  pllo[3] = 1.0;
250 
251  Char_t FunName[100];
252  sprintf(FunName, "FitfcnLG_%s", htoFit->GetName());
253 
254  langausFit = new TF1(FunName, HDQMUtil::langaufun, fr[0], fr[1], 4);
255  langausFit->SetParameters(sv);
256  langausFit->SetParNames("Width", "MP", "Area", "GSigma");
257 
258  for (int32_t i = 0; i < 4; i++) {
259  langausFit->SetParLimits(i, pllo[i], plhi[i]);
260  }
261 
262  try {
263  htoFit->Fit(langausFit, "R0"); // "R" fit in a range,"0" quiet fit
264 
265  langausFit->SetRange(fr[0], fr[1]);
266  langausFit->GetParameters(pLanGausS);
267  std::memcpy((void *)epLanGausS, (void *)langausFit->GetParErrors(), 4 * sizeof(double));
268 
269  chi2GausS = langausFit->GetChisquare(); // obtain chi^2
270  nDofGausS = langausFit->GetNDF(); // obtain ndf
271 
272  double sPeak, sFWHM;
273  HDQMUtil::langaupro(pLanGausS, sPeak, sFWHM);
274  pLanConv[0] = sPeak;
275  pLanConv[1] = sFWHM;
276  edm::LogInfo("fitUtility") << "langaupro: max " << sPeak << std::endl;
277  edm::LogInfo("fitUtility") << "langaupro: FWHM " << sFWHM << std::endl;
278  } catch (cms::Exception &iException) {
279  edm::LogError("fitUtility") << "problem in fitting " << htoFit->GetTitle()
280  << " \n\tDefault values of the parameters will be used";
281  pLanGausS[0] = -9999;
282  pLanGausS[1] = -9999;
283  pLanGausS[2] = -9999;
284  pLanGausS[3] = -9999;
285  epLanGausS[0] = -9999;
286  epLanGausS[1] = -9999;
287  epLanGausS[2] = -9999;
288  epLanGausS[3] = -9999;
289  pLanConv[0] = -9999;
290  pLanConv[1] = -9999;
291  chi2GausS = -9999;
292  nDofGausS = -9999;
293  }
294  } else {
295  pLanGausS[0] = -9999;
296  pLanGausS[1] = -9999;
297  pLanGausS[2] = -9999;
298  pLanGausS[3] = -9999;
299  epLanGausS[0] = -9999;
300  epLanGausS[1] = -9999;
301  epLanGausS[2] = -9999;
302  epLanGausS[3] = -9999;
303  pLanConv[0] = -9999;
304  pLanConv[1] = -9999;
305  chi2GausS = -9999;
306  nDofGausS = -9999;
307  }
308 
309  return htoFit->GetEntries();
310 }
int32_t langaupro(double *params, double &maxx, double &FWHM)
double langaufun(double *x, double *par)
double HDQMfitUtilities::getFitChi ( )
inline
int HDQMfitUtilities::getFitnDof ( )
inline
double HDQMfitUtilities::getGaussPar ( std::string  s)

Definition at line 428 of file HDQMfitUtilities.cc.

References pGausS.

Referenced by SiStripDQMHistoryHelper::setDBValuesForGauss().

428  {
429  if (s == "area")
430  return pGausS[0];
431  else if (s == "mean")
432  return pGausS[1];
433  else if (s == "sigma")
434  return pGausS[2];
435  else
436  return -99999;
437 }
double HDQMfitUtilities::getGaussParErr ( std::string  s)

Definition at line 439 of file HDQMfitUtilities.cc.

References epGausS.

439  {
440  if (s == "area")
441  return epGausS[0];
442  else if (s == "mean")
443  return epGausS[1];
444  else if (s == "sigma")
445  return epGausS[2];
446  else
447  return -99999;
448 }
double HDQMfitUtilities::getLanGaussConv ( std::string  s)

Definition at line 419 of file HDQMfitUtilities.cc.

References pLanConv.

Referenced by SiStripDQMHistoryHelper::setDBValuesForLandau().

419  {
420  if (s == "mpv")
421  return pLanConv[0];
422  else if (s == "fwhm")
423  return pLanConv[1];
424  else
425  return -99999;
426 }
double HDQMfitUtilities::getLanGaussPar ( std::string  s)

Definition at line 393 of file HDQMfitUtilities.cc.

References pLanGausS.

Referenced by SiStripDQMHistoryHelper::setDBValuesForLandau().

393  {
394  if (s == "landau_width")
395  return pLanGausS[0];
396  else if (s == "mpv")
397  return pLanGausS[1];
398  else if (s == "area")
399  return pLanGausS[2];
400  else if (s == "gauss_sigma")
401  return pLanGausS[3];
402  else
403  return -99999;
404 }
double HDQMfitUtilities::getLanGaussParErr ( std::string  s)

Definition at line 406 of file HDQMfitUtilities.cc.

References epLanGausS.

Referenced by SiStripDQMHistoryHelper::setDBValuesForLandau().

406  {
407  if (s == "landau_width")
408  return epLanGausS[0];
409  else if (s == "mpv")
410  return epLanGausS[1];
411  else if (s == "area")
412  return epLanGausS[2];
413  else if (s == "gauss_sigma")
414  return epLanGausS[3];
415  else
416  return -99999;
417 }
void HDQMfitUtilities::init ( void  )

Definition at line 176 of file HDQMfitUtilities.cc.

References chi2GausS, epGausS, epLanGausS, nDofGausS, pGausS, pLanConv, and pLanGausS.

Referenced by doGaussFit(), doLanGaussFit(), and HDQMfitUtilities().

176  {
177  pLanGausS[0] = 0;
178  pLanGausS[1] = 0;
179  pLanGausS[2] = 0;
180  pLanGausS[3] = 0;
181  epLanGausS[0] = 0;
182  epLanGausS[1] = 0;
183  epLanGausS[2] = 0;
184  epLanGausS[3] = 0;
185  pGausS[0] = 0;
186  pGausS[1] = 0;
187  pGausS[2] = 0;
188  epGausS[0] = 0;
189  epGausS[1] = 0;
190  epGausS[2] = 0;
191  pLanConv[0] = 0;
192  pLanConv[1] = 0;
193 
194  // if ( langausFit!=0 ) delete langausFit;
195  // if ( gausFit!=0 ) delete gausFit;
196 
197  chi2GausS = -9999.;
198  nDofGausS = -9999;
199 }

Member Data Documentation

double HDQMfitUtilities::chi2GausS
private

Definition at line 54 of file HDQMfitUtilities.h.

Referenced by doGaussFit(), doLanGaussFit(), and init().

double HDQMfitUtilities::epGausS[3]
private

Definition at line 52 of file HDQMfitUtilities.h.

Referenced by doGaussFit(), getGaussParErr(), and init().

double HDQMfitUtilities::epLanGausS[4]
private

Definition at line 51 of file HDQMfitUtilities.h.

Referenced by doLanGaussFit(), getLanGaussParErr(), and init().

TF1* HDQMfitUtilities::gausFit
private

Definition at line 57 of file HDQMfitUtilities.h.

Referenced by doGaussFit(), and ~HDQMfitUtilities().

TF1* HDQMfitUtilities::langausFit
private

Definition at line 56 of file HDQMfitUtilities.h.

Referenced by doGaussFit(), doLanGaussFit(), and ~HDQMfitUtilities().

int32_t HDQMfitUtilities::nDofGausS
private

Definition at line 55 of file HDQMfitUtilities.h.

Referenced by doGaussFit(), doLanGaussFit(), and init().

double HDQMfitUtilities::pGausS[3]
private

Definition at line 52 of file HDQMfitUtilities.h.

Referenced by doGaussFit(), getGaussPar(), and init().

double HDQMfitUtilities::pLanConv[2]
private

Definition at line 53 of file HDQMfitUtilities.h.

Referenced by doLanGaussFit(), getLanGaussConv(), and init().

double HDQMfitUtilities::pLanGausS[4]
private

Definition at line 51 of file HDQMfitUtilities.h.

Referenced by doLanGaussFit(), getLanGaussPar(), and init().