CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Attributes
DTMeanTimerFitter Class Reference

#include <DTMeanTimerFitter.h>

Public Member Functions

 DTMeanTimerFitter (TFile *file)
 Constructor. More...
 
std::vector< float > evaluateVDriftAndReso (const TString &N)
 Fit the TMax histos and evaluate VDrift and resolution. More...
 
TF1 * fitTMax (TH1F *histo)
 Really do the fit. More...
 
void setVerbosity (unsigned int lvl)
 Set the verbosity of the output: 0 = silent, 1 = info, 2 = debug. More...
 
virtual ~DTMeanTimerFitter ()
 Destructor. More...
 

Private Attributes

TFile * hDebugFile
 
TFile * hInputFile
 
unsigned int theVerbosityLevel
 

Detailed Description

Fit the Tmax histograms with a gaussian returning the mean values and the sigmas.

Author
S. Bolognesi - INFN Torino

Definition at line 17 of file DTMeanTimerFitter.h.

Constructor & Destructor Documentation

DTMeanTimerFitter::DTMeanTimerFitter ( TFile *  file)

Constructor.

Definition at line 21 of file DTMeanTimerFitter.cc.

References hDebugFile.

22  //hInputFile = new TFile("DTTMaxHistos.root", "READ");
23  hDebugFile = new TFile("DTMeanTimerFitter.root", "RECREATE");
24 }
unsigned int theVerbosityLevel
DTMeanTimerFitter::~DTMeanTimerFitter ( )
virtual

Destructor.

Definition at line 26 of file DTMeanTimerFitter.cc.

References hDebugFile.

26 { hDebugFile->Close(); }

Member Function Documentation

vector< float > DTMeanTimerFitter::evaluateVDriftAndReso ( const TString &  N)

Fit the TMax histos and evaluate VDrift and resolution.

Definition at line 28 of file DTMeanTimerFitter.cc.

References submitPVResolutionJobs::count, counter, cppFunctionSkipper::exception, validate-o2o-wbm::f1, fitTMax(), hTMaxCell::GetT0Factor(), hDebugFile, hInputFile, mergeVDriftHistosByStation::histos, hTMaxCell::hTmax123, hTMaxCell::hTmax124s72, hTMaxCell::hTmax124s78, hTMaxCell::hTmax134s72, hTMaxCell::hTmax134s78, hTMaxCell::hTmax234, hTMaxCell::hTmax_0, hTMaxCell::hTmax_2t0, hTMaxCell::hTmax_3t0, hTMaxCell::hTmax_t0, mps_fire::i, SiStripPI::mean, and mathSSE::sqrt().

28  {
29  // Retrieve histogram sets
31  vector<float> vDriftAndReso;
32 
33  // Check that the histo for this cell exists
34  if (histos->hTmax123 != nullptr) {
35  vector<TH1F *> hTMax; // histograms for <T_max> calculation
36  vector<TH1F *> hT0; // histograms for T0 evaluation
37  hTMax.push_back(histos->hTmax123);
38  hTMax.push_back(histos->hTmax124s72);
39  hTMax.push_back(histos->hTmax124s78);
40  hTMax.push_back(histos->hTmax134s72);
41  hTMax.push_back(histos->hTmax134s78);
42  hTMax.push_back(histos->hTmax234);
43 
44  hT0.push_back(histos->hTmax_3t0);
45  hT0.push_back(histos->hTmax_2t0);
46  hT0.push_back(histos->hTmax_t0);
47  hT0.push_back(histos->hTmax_0);
48 
49  vector<Double_t> factor; // factor relating the width of the Tmax distribution
50  // and the cell resolution
51  factor.push_back(sqrt(2. / 3.)); // hTmax123
52  factor.push_back(sqrt(2. / 7.)); // hTmax124s72
53  factor.push_back(sqrt(8. / 7.)); // hTmax124s78
54  factor.push_back(sqrt(2. / 7.)); // hTmax134s72
55  factor.push_back(sqrt(8. / 7.)); // hTmax134s78
56  factor.push_back(sqrt(2. / 3.)); // hTmax234
57 
58  // Retrieve the gaussian mean and sigma for each TMax histogram
59  vector<Double_t> mean;
60  vector<Double_t> sigma;
61  vector<Double_t> count; //number of entries
62 
63  for (vector<TH1F *>::const_iterator ith = hTMax.begin(); ith != hTMax.end(); ++ith) {
64  TF1 *funct = fitTMax(*ith);
65  if (!funct) {
66  edm::LogError("DTMeanTimerFitter") << "Error when fitting TMax..histogram name" << (*ith)->GetName();
67  // return empty or -1 filled vector?
68  vector<float> defvec(6, -1);
69  return defvec;
70  }
71  hDebugFile->cd();
72  (*ith)->Write();
73 
74  // Get mean, sigma and number of entries of each histogram
75  mean.push_back(funct->GetParameter(1));
76  sigma.push_back(funct->GetParameter(2));
77  count.push_back((*ith)->GetEntries());
78  }
79 
80  Double_t tMaxMean = 0.;
81  Double_t wTMaxSum = 0.;
82  Double_t sigmaT = 0.;
83  Double_t wSigmaSum = 0.;
84 
85  //calculate total mean and sigma
86  for (int i = 0; i <= 5; i++) {
87  if (count[i] < 200)
88  continue;
89  tMaxMean += mean[i] * (count[i] / (sigma[i] * sigma[i]));
90  wTMaxSum += count[i] / (sigma[i] * sigma[i]);
91  sigmaT += count[i] * factor[i] * sigma[i];
92  wSigmaSum += count[i];
93  // cout << "TMaxMean "<<i<<": "<< mean[i] << " entries: " << count[i]
94  // << " sigma: " << sigma[i]
95  // << " weight: " << (count[i]/(sigma[i]*sigma[i])) << endl;
96  }
97  if ((!wTMaxSum) || (!wSigmaSum)) {
98  edm::LogError("DTMeanTimerFitter") << "Error zero sum of weights..returning default";
99  vector<float> defvec(6, -1);
100  return defvec;
101  }
102 
103  tMaxMean /= wTMaxSum;
104  sigmaT /= wSigmaSum;
105 
106  //calculate v_drift and resolution
107  Double_t vDrift = 2.1 / tMaxMean; //2.1 is the half cell length in cm
108  Double_t reso = vDrift * sigmaT;
109  vDriftAndReso.push_back(vDrift);
110  vDriftAndReso.push_back(reso);
111  //if(theVerbosityLevel >= 1)
112  edm::LogVerbatim("DTMeanTimerFitter")
113  << " final TMaxMean=" << tMaxMean << " sigma= " << sigmaT << " v_d and reso: " << vDrift << " " << reso << endl;
114 
115  // Order t0 histogram by number of entries (choose histograms with higher nr. of entries)
116  map<Double_t, TH1F *> hEntries;
117  for (vector<TH1F *>::const_iterator ith = hT0.begin(); ith != hT0.end(); ++ith) {
118  hEntries[(*ith)->GetEntries()] = (*ith);
119  }
120 
121  // add at the end of hT0 the two hists with the higher number of entries
122  int counter = 0;
123  for (map<Double_t, TH1F *>::reverse_iterator iter = hEntries.rbegin(); iter != hEntries.rend(); ++iter) {
124  counter++;
125  if (counter == 1)
126  hT0.push_back(iter->second);
127  else if (counter == 2) {
128  hT0.push_back(iter->second);
129  break;
130  }
131  }
132 
133  // Retrieve the gaussian mean and sigma of histograms for Delta(t0) evaluation
134  vector<Double_t> meanT0;
135  vector<Double_t> sigmaT0;
136  vector<Double_t> countT0;
137 
138  for (vector<TH1F *>::const_iterator ith = hT0.begin(); ith != hT0.end(); ++ith) {
139  try {
140  (*ith)->Fit("gaus");
141  } catch (std::exception const &) {
142  edm::LogError("DTMeanTimerFitter") << "Exception when fitting T0..histogram " << (*ith)->GetName();
143  // return empty or -1 filled vector?
144  vector<float> defvec(6, -1);
145  return defvec;
146  }
147  TF1 *f1 = (*ith)->GetFunction("gaus");
148  // Get mean, sigma and number of entries of the histograms
149  meanT0.push_back(f1->GetParameter(1));
150  sigmaT0.push_back(f1->GetParameter(2));
151  countT0.push_back((*ith)->GetEntries());
152  }
153  //calculate Delta(t0)
154  if (hT0.size() != 6) { // check if you have all the t0 hists
155  edm::LogVerbatim("DTMeanTimerFitter") << "t0 histograms = " << hT0.size();
156  for (int i = 1; i <= 4; i++) {
157  vDriftAndReso.push_back(-1);
158  }
159  return vDriftAndReso;
160  }
161 
162  for (int it0 = 0; it0 <= 2; it0++) {
163  if ((countT0[it0] > 200) && (countT0[it0 + 1] > 200)) {
164  Double_t deltaT0 = meanT0[it0] - meanT0[it0 + 1];
165  vDriftAndReso.push_back(deltaT0);
166  } else
167  vDriftAndReso.push_back(999.);
168  }
169  //deltat0 using hists with max nr. of entries
170  if ((countT0[4] > 200) && (countT0[5] > 200)) {
171  Double_t t0Diff = histos->GetT0Factor(hT0[4]) - histos->GetT0Factor(hT0[5]);
172  Double_t deltaT0MaxEntries = (meanT0[4] - meanT0[5]) / t0Diff;
173  vDriftAndReso.push_back(deltaT0MaxEntries);
174  } else
175  vDriftAndReso.push_back(999.);
176  } else {
177  for (int i = 1; i <= 6; i++) {
178  // 0=vdrift, 1=reso, 2=(3deltat0-2deltat0), 3=(2deltat0-1deltat0),
179  // 4=(1deltat0-0deltat0), 5=deltat0 from hists with max entries,
180  vDriftAndReso.push_back(-1);
181  }
182  }
183  return vDriftAndReso;
184 }
Log< level::Info, true > LogVerbatim
TH1F * hTmax134s78
Definition: vDriftHistos.h:514
TH1F * hTmax_2t0
Definition: vDriftHistos.h:523
TH1F * hTmax_3t0
Definition: vDriftHistos.h:516
int GetT0Factor(TH1F *hist)
Definition: vDriftHistos.h:495
TH1F * hTmax234
Definition: vDriftHistos.h:515
Log< level::Error, false > LogError
TH1F * hTmax_t0
Definition: vDriftHistos.h:530
TH1F * hTmax134s72
Definition: vDriftHistos.h:513
TH1F * hTmax124s72
Definition: vDriftHistos.h:511
TF1 * fitTMax(TH1F *histo)
Really do the fit.
TH1F * hTmax_0
Definition: vDriftHistos.h:537
T sqrt(T t)
Definition: SSEVec.h:19
TH1F * hTmax124s78
Definition: vDriftHistos.h:512
#define N
Definition: blowfish.cc:9
static std::atomic< unsigned int > counter
TH1F * hTmax123
Definition: vDriftHistos.h:510
TF1 * DTMeanTimerFitter::fitTMax ( TH1F *  histo)

Really do the fit.

Definition at line 186 of file DTMeanTimerFitter.cc.

References cppFunctionSkipper::exception, validate-o2o-wbm::f1, LogDebug, and sistrip::SpyUtilities::range().

Referenced by evaluateVDriftAndReso().

186  {
187  // Find distribution peak and fit range
188  Double_t peak = (((((histo->GetXaxis())->GetXmax()) - ((histo->GetXaxis())->GetXmin())) / histo->GetNbinsX()) *
189  (histo->GetMaximumBin())) +
190  ((histo->GetXaxis())->GetXmin());
191  //if(theVerbosityLevel >= 1)
192  LogDebug("DTMeanTimerFitter") << "Peak " << peak << " : "
193  << "xmax " << ((histo->GetXaxis())->GetXmax()) << " xmin "
194  << ((histo->GetXaxis())->GetXmin()) << " nbin " << histo->GetNbinsX()
195  << " bin with max " << (histo->GetMaximumBin());
196  Double_t range = 2. * histo->GetRMS();
197 
198  // Fit each Tmax histogram with a Gaussian in a restricted interval
199  TF1 *rGaus = new TF1("rGaus", "gaus", peak - range, peak + range);
200  rGaus->SetMarkerSize(); // just silence gcc complainining about unused var
201  try {
202  histo->Fit("rGaus", "R");
203  } catch (std::exception const &) {
204  edm::LogError("DTMeanTimerFitter") << "Exception when fitting TMax..histogram " << histo->GetName()
205  << " setting return function pointer to zero";
206  return nullptr;
207  }
208  TF1 *f1 = histo->GetFunction("rGaus");
209  return f1;
210 }
Log< level::Error, false > LogError
const uint16_t range(const Frame &aFrame)
#define LogDebug(id)
void DTMeanTimerFitter::setVerbosity ( unsigned int  lvl)
inline

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

Definition at line 29 of file DTMeanTimerFitter.h.

References theVerbosityLevel.

29 { theVerbosityLevel = lvl; }
unsigned int theVerbosityLevel

Member Data Documentation

TFile* DTMeanTimerFitter::hDebugFile
private
TFile* DTMeanTimerFitter::hInputFile
private

Definition at line 37 of file DTMeanTimerFitter.h.

Referenced by evaluateVDriftAndReso().

unsigned int DTMeanTimerFitter::theVerbosityLevel
private

Definition at line 39 of file DTMeanTimerFitter.h.

Referenced by setVerbosity().