CMS 3D CMS Logo

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

#include <DTTimeBoxFitter.h>

Public Member Functions

 DTTimeBoxFitter (const TString &debugFileName=TString(""))
 Constructor. More...
 
std::pair< double, double > fitTimeBox (TH1F *hTimeBox)
 Fit the rising edge of the time box returning mean value and sigma (first and second respectively) More...
 
void getFitSeeds (TH1F *hTBox, double &mean, double &sigma, double &tBoxMax, double &xFitMin, double &xFitMax)
 Automatically compute the seeds the range to be used for time box fit. More...
 
void getInteractiveFitSeeds (TH1F *hTBox, double &mean, double &sigma, double &tBoxMax, double &xFitMin, double &xFitMax)
 Ask the user to provide the seeds. More...
 
void setFitSigma (double sigma)
 
void setInteractiveFit (bool isInteractive)
 Switch to interactive fit. More...
 
void setRebinning (int reb)
 Set the rebin. More...
 
void setVerbosity (unsigned int lvl)
 Set the verbosity of the output: 0 = silent, 1 = info, 2 = debug. More...
 
virtual ~DTTimeBoxFitter ()
 Destructor. More...
 

Private Attributes

TFile * hDebugFile
 
bool interactiveFit
 
int rebin
 
double theSigma
 
unsigned int theVerbosityLevel
 

Detailed Description

Fit the rising edge of the time box with the integral of a gaussian returning the mean value and the sigma.

Author
G. Cerminara - INFN Torino

Definition at line 17 of file DTTimeBoxFitter.h.

Constructor & Destructor Documentation

◆ DTTimeBoxFitter()

DTTimeBoxFitter::DTTimeBoxFitter ( const TString &  debugFileName = TString(""))

Constructor.

Definition at line 21 of file DTTimeBoxFitter.cc.

References hDebugFile, interactiveFit, and rebin.

22  : hDebugFile(nullptr), theVerbosityLevel(0), theSigma(10.) {
23  // Create a root file for debug output only if needed
24  if (debugFileName != "")
25  hDebugFile = new TFile(debugFileName.Data(), "RECREATE");
26  interactiveFit = false;
27  rebin = 1;
28 }
unsigned int theVerbosityLevel

◆ ~DTTimeBoxFitter()

DTTimeBoxFitter::~DTTimeBoxFitter ( )
virtual

Destructor.

Definition at line 30 of file DTTimeBoxFitter.cc.

References hDebugFile.

30  {
31  if (hDebugFile != nullptr)
32  hDebugFile->Close();
33 }

Member Function Documentation

◆ fitTimeBox()

pair< double, double > DTTimeBoxFitter::fitTimeBox ( TH1F *  hTimeBox)

Fit the rising edge of the time box returning mean value and sigma (first and second respectively)

Compute the ttrig (in ns) from the Time Box.

Definition at line 36 of file DTTimeBoxFitter.cc.

References simKBmtfDigis_cfi::chiSquare, gather_cfg::cout, getFitSeeds(), getInteractiveFitSeeds(), interactiveFit, intGauss(), SiStripPI::mean, fileinputsource_cfi::option, rebin, and theVerbosityLevel.

Referenced by DTTTrigWriter::analyze().

36  {
37  hTimeBox->Rebin(rebin);
38 
39  // Check if the histo contains any entry
40  if (hTimeBox->GetEntries() == 0) {
41  cout << "[DTTimeBoxFitter]***Error: the time box contains no entry!" << endl;
42  return make_pair(-1, -1);
43  }
44 
45  if (hTimeBox->GetEntries() < 50000) {
46  hTimeBox->Rebin(2);
47  }
48 
49  // Get seeds for the fit
50  // The TimeBox range to be fitted (the rising edge)
51  double xFitMin = 0; // Min value for the fit
52  double xFitMax = 0; // Max value for the fit
53  double xValue = 0; // The mean value of the gaussian
54  double xFitSigma = 0; // The sigma of the gaussian
55  double tBoxMax = 0; // The max of the time box, it is used as seed for gaussian integral
56 
57  //hTimeBox->Rebin(2); //FIXME: Temporary for low statistics
58 
59  TH1F* hTimeBoxForSeed = (TH1F*)hTimeBox->Clone(); //FIXME: test
60  if (!interactiveFit) {
61  getFitSeeds(hTimeBoxForSeed, xValue, xFitSigma, tBoxMax, xFitMin, xFitMax);
62  } else {
63  getInteractiveFitSeeds(hTimeBoxForSeed, xValue, xFitSigma, tBoxMax, xFitMin, xFitMax);
64  }
65 
66  // Define the fitting function and use fit seeds
67  TF1* fIntGaus = new TF1("IntGauss", intGauss, xFitMin, xFitMax, 3);
68  fIntGaus->SetParName(0, "Constant");
69  fIntGaus->SetParameter(0, tBoxMax);
70  fIntGaus->SetParName(1, "Mean");
71  fIntGaus->SetParameter(1, xValue);
72  fIntGaus->SetParName(2, "Sigma");
73  fIntGaus->SetParameter(2, xFitSigma);
74  fIntGaus->SetLineColor(kRed);
75 
76  // Fit the histo
77  string option = "Q";
78  if (theVerbosityLevel >= 2)
79  option = "";
80 
81  hTimeBox->Fit(fIntGaus, option.c_str(), "", xFitMin, xFitMax);
82 
83  // Get fitted parameters
84  double mean = fIntGaus->GetParameter("Mean");
85  double sigma = fIntGaus->GetParameter("Sigma");
86  // double constant = fIntGaus->GetParameter("Constant");
87  double chiSquare = fIntGaus->GetChisquare() / fIntGaus->GetNDF();
88 
89  if (theVerbosityLevel >= 1) {
90  cout << " === Fit Results: " << endl;
91  cout << " Fitted mean = " << mean << endl;
92  cout << " Fitted sigma = " << sigma << endl;
93  cout << " Reduced Chi Square = " << chiSquare << endl;
94  }
95  return make_pair(mean, sigma);
96 }
double intGauss(double *x, double *par)
void getFitSeeds(TH1F *hTBox, double &mean, double &sigma, double &tBoxMax, double &xFitMin, double &xFitMax)
Automatically compute the seeds the range to be used for time box fit.
void getInteractiveFitSeeds(TH1F *hTBox, double &mean, double &sigma, double &tBoxMax, double &xFitMin, double &xFitMax)
Ask the user to provide the seeds.
unsigned int theVerbosityLevel

◆ getFitSeeds()

void DTTimeBoxFitter::getFitSeeds ( TH1F *  hTBox,
double &  mean,
double &  sigma,
double &  tBoxMax,
double &  xFitMin,
double &  xFitMax 
)

Automatically compute the seeds the range to be used for time box fit.

Definition at line 125 of file DTTimeBoxFitter.cc.

References funct::abs(), filterCSVwithJSON::copy, gather_cfg::cout, dumpMFGeometry_cfg::delta, change_name::diff, HLT_2024v14_cff::gap, hDebugFile, mps_fire::i, createfilelist::int, MergeJob_cfg::interval, dqmiolumiharvest::j, SiStripPI::mean, ecalPedestalPCLworker_cfi::nBins, GetRecoTauVFromDQM_MC_cff::next, theSigma, theVerbosityLevel, DiMuonV_cfg::threshold, trackerHitRTTI::vector, multiplicitycorr_cfi::xMax, and photonAnalyzer_cfi::xMin.

Referenced by fitTimeBox().

126  {
127  if (theVerbosityLevel >= 1)
128  cout << " === Looking for fit seeds in Time Box:" << endl;
129 
130  // The approximate width of the time box
131  static const int tBoxWidth = 400; //FIXE: tune it
132 
133  int nBins = hTBox->GetNbinsX();
134  const int xMin = (int)hTBox->GetXaxis()->GetXmin();
135  const int xMax = (int)hTBox->GetXaxis()->GetXmax();
136  const int nEntries = (int)hTBox->GetEntries();
137 
138  double binValue = (double)(xMax - xMin) / (double)nBins;
139 
140  // Compute a threshold for TimeBox discrimination
141  const double threshold = binValue * nEntries / (double)(tBoxWidth * 3.);
142  if (theVerbosityLevel >= 2)
143  cout << " Threshold for logic time box is (# entries): " << threshold << endl;
144 
145  int nRebins = 0; // protection for infinite loop
146  while (threshold > hTBox->GetMaximum() / 2. && nRebins < 5) {
147  cout << " Rebinning!" << endl;
148  hTBox->Rebin(2);
149  nBins = hTBox->GetNbinsX();
150  binValue = (double)(xMax - xMin) / (double)nBins;
151  nRebins++;
152  }
153 
154  if (hDebugFile != nullptr)
155  hDebugFile->cd();
156  TString hLName = TString(hTBox->GetName()) + "L";
157  TH1F hLTB(hLName.Data(), "Logic Time Box", nBins, xMin, xMax);
158  // Loop over all time box bins and discriminate them accordigly to the threshold
159  for (int i = 1; i <= nBins; i++) {
160  if (hTBox->GetBinContent(i) > threshold)
161  hLTB.SetBinContent(i, 1);
162  }
163  if (hDebugFile != nullptr)
164  hLTB.Write();
165 
166  // Look for the time box in the "logic histo" and save beginning and lenght of each plateau
167  vector<pair<int, int> > startAndLenght;
168  if (theVerbosityLevel >= 2)
169  cout << " Look for rising and folling edges of logic time box: " << endl;
170  int start = -1;
171  int lenght = -1;
172  for (int j = 1; j < nBins; j++) {
173  int diff = (int)hLTB.GetBinContent(j + 1) - (int)hLTB.GetBinContent(j);
174  if (diff == 1) { // This is a rising edge
175  start = j;
176  lenght = 1;
177  if (theVerbosityLevel >= 2) {
178  cout << " >>>----" << endl;
179  cout << " bin: " << j << " is a rising edge" << endl;
180  }
181  } else if (diff == -1) { // This is a falling edge
182  startAndLenght.push_back(make_pair(start, lenght));
183  if (theVerbosityLevel >= 2) {
184  cout << " bin: " << j << " is a falling edge, lenght is: " << lenght << endl;
185  cout << " <<<----" << endl;
186  }
187  start = -1;
188  lenght = -1;
189  } else if (diff == 0 && (int)hLTB.GetBinContent(j) == 1) { // This is a bin within the plateau
190  lenght++;
191  }
192  }
193 
194  if (theVerbosityLevel >= 2)
195  cout << " Add macro intervals: " << endl;
196  // Look for consecutive plateau: 2 plateau are consecutive if the gap is < 5 ns
197  vector<pair<int, int> > superIntervals;
198  for (vector<pair<int, int> >::const_iterator interval = startAndLenght.begin(); interval != startAndLenght.end();
199  ++interval) {
200  pair<int, int> theInterval = (*interval);
201  vector<pair<int, int> >::const_iterator next = interval;
202  while (++next != startAndLenght.end()) {
203  int gap = (*next).first - (theInterval.first + theInterval.second);
204  double gabInNs = binValue * gap;
205  if (theVerbosityLevel >= 2)
206  cout << " gap: " << gabInNs << "(ns)" << endl;
207  if (gabInNs > 20) {
208  break;
209  } else {
210  theInterval = make_pair(theInterval.first, theInterval.second + gap + (*next).second);
211  superIntervals.push_back(theInterval);
212  if (theVerbosityLevel >= 2)
213  cout << " Add interval, start: " << theInterval.first << " width: " << theInterval.second << endl;
214  }
215  }
216  }
217 
218  // merge the vectors of intervals
219  copy(superIntervals.begin(), superIntervals.end(), back_inserter(startAndLenght));
220 
221  // Look for the plateau of the right lenght
222  if (theVerbosityLevel >= 2)
223  cout << " Look for the best interval:" << endl;
224  int delta = 999999;
225  int beginning = -1;
226  int tbWidth = -1;
227  for (vector<pair<int, int> >::const_iterator stAndL = startAndLenght.begin(); stAndL != startAndLenght.end();
228  ++stAndL) {
229  if (abs((*stAndL).second - tBoxWidth) < delta) {
230  delta = abs((*stAndL).second - tBoxWidth);
231  beginning = (*stAndL).first;
232  tbWidth = (*stAndL).second;
233  if (theVerbosityLevel >= 2)
234  cout << " Candidate: First: " << beginning << ", width: " << tbWidth << ", delta: " << delta << endl;
235  }
236  }
237 
238  mean = xMin + beginning * binValue;
239  sigma = theSigma; //FIXME: estimate it!
240 
241  tBoxMax = hTBox->GetMaximum();
242 
243  // Define the fit range
244  xFitMin = mean - 5. * sigma;
245  xFitMax = mean + 5. * sigma;
246 
247  if (theVerbosityLevel >= 1) {
248  cout << " Time Box Rising edge: " << mean << endl;
249  cout << " Time Box Width: " << tbWidth * binValue << endl;
250  cout << " = Seeds and range for fit:" << endl;
251  cout << " Seed mean = " << mean << endl;
252  cout << " Seed sigma = " << sigma << endl;
253  cout << " Fitting from = " << xFitMin << " to " << xFitMax << endl << endl;
254  }
255 }
Definition: start.py:1
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int theVerbosityLevel

◆ getInteractiveFitSeeds()

void DTTimeBoxFitter::getInteractiveFitSeeds ( TH1F *  hTBox,
double &  mean,
double &  sigma,
double &  tBoxMax,
double &  xFitMin,
double &  xFitMax 
)

Ask the user to provide the seeds.

Definition at line 99 of file DTTimeBoxFitter.cc.

References gather_cfg::cout, SiStripPI::mean, theSigma, and theVerbosityLevel.

Referenced by fitTimeBox().

100  {
101  if (theVerbosityLevel >= 1)
102  cout << " === Insert seeds for the Time Box fit:" << endl;
103 
104  cout << "Inser the fit mean:" << endl;
105  cin >> mean;
106 
107  sigma = theSigma; //FIXME: estimate it!
108 
109  tBoxMax = hTBox->GetMaximum();
110 
111  // Define the fit range
112  xFitMin = mean - 5. * sigma;
113  xFitMax = mean + 5. * sigma;
114 
115  if (theVerbosityLevel >= 1) {
116  cout << " Time Box Rising edge: " << mean << endl;
117  cout << " = Seeds and range for fit:" << endl;
118  cout << " Seed mean = " << mean << endl;
119  cout << " Seed sigma = " << sigma << endl;
120  cout << " Fitting from = " << xFitMin << " to " << xFitMax << endl << endl;
121  }
122 }
unsigned int theVerbosityLevel

◆ setFitSigma()

void DTTimeBoxFitter::setFitSigma ( double  sigma)
inline

Definition at line 46 of file DTTimeBoxFitter.h.

References theSigma.

Referenced by DTTTrigWriter::DTTTrigWriter().

46 { theSigma = sigma; }

◆ setInteractiveFit()

void DTTimeBoxFitter::setInteractiveFit ( bool  isInteractive)
inline

Switch to interactive fit.

Definition at line 41 of file DTTimeBoxFitter.h.

References interactiveFit.

41 { interactiveFit = isInteractive; }

◆ setRebinning()

void DTTimeBoxFitter::setRebinning ( int  reb)
inline

Set the rebin.

Definition at line 44 of file DTTimeBoxFitter.h.

References rebin.

44 { rebin = reb; }

◆ setVerbosity()

void DTTimeBoxFitter::setVerbosity ( unsigned int  lvl)
inline

Set the verbosity of the output: 0 = silent, 1 = info, 2 = debug.

Definition at line 38 of file DTTimeBoxFitter.h.

References theVerbosityLevel.

Referenced by DTTTrigWriter::DTTTrigWriter().

38 { theVerbosityLevel = lvl; }
unsigned int theVerbosityLevel

Member Data Documentation

◆ hDebugFile

TFile* DTTimeBoxFitter::hDebugFile
private

Definition at line 50 of file DTTimeBoxFitter.h.

Referenced by DTTimeBoxFitter(), getFitSeeds(), and ~DTTimeBoxFitter().

◆ interactiveFit

bool DTTimeBoxFitter::interactiveFit
private

Definition at line 53 of file DTTimeBoxFitter.h.

Referenced by DTTimeBoxFitter(), fitTimeBox(), and setInteractiveFit().

◆ rebin

int DTTimeBoxFitter::rebin
private

Definition at line 54 of file DTTimeBoxFitter.h.

Referenced by DTTimeBoxFitter(), fitTimeBox(), and setRebinning().

◆ theSigma

double DTTimeBoxFitter::theSigma
private

Definition at line 55 of file DTTimeBoxFitter.h.

Referenced by getFitSeeds(), getInteractiveFitSeeds(), and setFitSigma().

◆ theVerbosityLevel

unsigned int DTTimeBoxFitter::theVerbosityLevel
private

Definition at line 52 of file DTTimeBoxFitter.h.

Referenced by fitTimeBox(), getFitSeeds(), getInteractiveFitSeeds(), and setVerbosity().