CMS 3D CMS Logo

LA_Fitter.cc
Go to the documentation of this file.
3 
4 #include <cmath>
5 #include <boost/algorithm/string/erase.hpp>
6 #include <TF1.h>
7 
9  for (Book::iterator it = book.begin(".*" + method(WIDTH)); it != book.end(); ++it) {
10  it->second->SetTitle("Mean Cluster Width;tan#theta_{t}");
11  TH1* const p = it->second;
12  if (p->GetEntries() < 400) {
13  delete p;
14  book[it->first] = nullptr;
15  continue;
16  }
17  p->SetTitle(";tan#theta_{t};");
18  const float min = p->GetMinimum();
19  const float max = p->GetMaximum();
20  float xofmin = p->GetBinCenter(p->GetMinimumBin());
21  if (xofmin > 0.0 || xofmin < -0.15)
22  xofmin = -0.05;
23  const float xofmax = p->GetBinCenter(p->GetMaximumBin());
24 
25  TF1* const fit = new TF1("LA_profile_fit", "[2]*(TMath::Abs(x-[0]))+[1]", -1, 1);
26  fit->SetParLimits(0, -0.15, 0.01);
27  fit->SetParLimits(1, 0.6 * min, 1.25 * min);
28  fit->SetParLimits(2, 0.1, 10);
29  fit->SetParameters(xofmin, min, (max - min) / fabs(xofmax - xofmin));
30 
31  int badfit = p->Fit(fit, "IEQ", "", -.5, .3);
32  if (badfit)
33  badfit = p->Fit(fit, "IEQ", "", -.46, .26);
34  if (badfit) {
35  book.erase(it);
36  }
37  }
38 }
39 
41  for (Book::iterator it = book.begin(".*_all"); it != book.end(); ++it) {
42  const std::string base = boost::erase_all_copy(it->first, "_all");
43 
44  std::vector<Book::iterator> rebin_hists;
45  const Book::iterator& all = it;
46  rebin_hists.push_back(all);
47  Book::iterator w1 = book.find(base + "_w1");
48  rebin_hists.push_back(w1);
49  Book::iterator var_w2 = book.find(base + method(AVGV2, false));
50  rebin_hists.push_back(var_w2);
51  Book::iterator var_w3 = book.find(base + method(AVGV3, false));
52  rebin_hists.push_back(var_w3);
53 
54  const unsigned rebin = std::max(var_w2 == book.end() ? 0 : find_rebin(var_w2->second),
55  var_w3 == book.end() ? 0 : find_rebin(var_w3->second));
56  for (const auto& it : rebin_hists)
57  if (it != book.end())
58  it->second->Rebin(rebin > 1 ? rebin < 7 ? rebin : 6 : 1);
59 
60  TH1* const prob_w1 =
61  w1 == book.end() ? nullptr : subset_probability(base + method(PROB1, false), w1->second, all->second);
62  TH1* const rmsv_w2 =
63  var_w2 == book.end() ? nullptr : rms_profile(base + method(RMSV2, false), (TProfile* const)var_w2->second);
64  TH1* const rmsv_w3 =
65  var_w3 == book.end() ? nullptr : rms_profile(base + method(RMSV3, false), (TProfile* const)var_w3->second);
66 
67  std::vector<TH1*> fit_hists;
68  if (prob_w1) {
69  book.book(base + method(PROB1, false), prob_w1);
70  fit_hists.push_back(prob_w1);
71  prob_w1->SetTitle("Width==1 Probability;tan#theta_{t}-(dx/dz)_{reco}");
72  }
73  if (var_w2 != book.end()) {
74  book.book(base + method(RMSV2, false), rmsv_w2);
75  fit_hists.push_back(var_w2->second);
76  var_w2->second->SetTitle("Width==2 Mean Variance;tan#theta_{t}-(dx/dz)_{reco}");
77  fit_hists.push_back(rmsv_w2);
78  rmsv_w2->SetTitle("Width==2 RMS Variance;tan#theta_{t}-(dx/dz)_{reco}");
79  }
80  if (var_w3 != book.end()) {
81  book.book(base + method(RMSV3, false), rmsv_w3);
82  fit_hists.push_back(var_w3->second);
83  var_w3->second->SetTitle("Width==3 Mean Variance;tan#theta_{t}-(dx/dz)_{reco}");
84  fit_hists.push_back(rmsv_w3);
85  rmsv_w3->SetTitle("Width==3 RMS Variance;tan#theta_{t}-(dx/dz)_{reco}");
86  }
87 
88  if (fit_hists.empty())
89  continue;
90  const unsigned bins = fit_hists[0]->GetNbinsX();
91  const unsigned guess = fit_hists[0]->FindBin(0);
92  const std::pair<unsigned, unsigned> range(guess - bins / 30, guess + bins / 30);
93 
94  for (auto const& hist : fit_hists) {
96  if (chi2) {
97  book.book(chi2->GetName(), chi2);
98  chi2->SetTitle("Symmetry #chi^{2};tan#theta_{t}-(dx/dz)_{reco}");
99  }
100  }
101  }
102 }
103 
104 unsigned LA_Filler_Fitter::find_rebin(const TH1* const hist) {
105  const double mean = hist->GetMean();
106  const double rms = hist->GetRMS();
107  const int begin = std::min(1, hist->GetXaxis()->FindFixBin(mean - rms));
108  const int end = std::max(hist->GetNbinsX(), hist->GetXaxis()->FindFixBin(mean + rms)) + 1;
109  unsigned current_hole(0), max_hole(0);
110  for (int i = begin; i < end; i++) {
111  if (!hist->GetBinError(i))
112  current_hole++;
113  else if (current_hole) {
114  max_hole = std::max(current_hole, max_hole);
115  current_hole = 0;
116  }
117  }
118  return max_hole + 1;
119 }
120 
121 TH1* LA_Filler_Fitter::rms_profile(const std::string name, const TProfile* const prof) {
122  const int bins = prof->GetNbinsX();
123  TH1* const rms =
124  new TH1F(name.c_str(), "", bins, prof->GetBinLowEdge(1), prof->GetBinLowEdge(bins) + prof->GetBinWidth(bins));
125  for (int i = 1; i <= bins; i++) {
126  const double Me = prof->GetBinError(i);
127  const double neff = prof->GetBinEntries(
128  i); //Should be prof->GetBinEffectiveEntries(i);, not availible this version ROOT. This is only ok for unweighted fills
129  rms->SetBinContent(i, Me * sqrt(neff));
130  rms->SetBinError(i, Me / sqrt(2));
131  }
132  return rms;
133 }
134 
135 TH1* LA_Filler_Fitter::subset_probability(const std::string name, const TH1* const subset, const TH1* const total) {
136  const int bins = subset->GetNbinsX();
137  TH1* const prob = new TH1F(
138  name.c_str(), "", bins, subset->GetBinLowEdge(1), subset->GetBinLowEdge(bins) + subset->GetBinWidth(bins));
139  for (int i = 1; i <= bins; i++) {
140  const double s = subset->GetBinContent(i);
141  const double T = total->GetBinContent(i);
142  const double B = T - s;
143 
144  const double p = T ? s / T : 0;
145  const double perr = T ? ((s && B) ? sqrt(s * s * B + B * B * s) / (T * T) : 1 / T) : 0;
146 
147  prob->SetBinContent(i, p);
148  prob->SetBinError(i, perr);
149  }
150  return prob;
151 }
FastTimerService_cff.range
range
Definition: FastTimerService_cff.py:34
mps_fire.i
i
Definition: mps_fire.py:428
SiStripPI::mean
Definition: SiStripPayloadInspectorHelper.h:169
LA_Filler_Fitter::find_rebin
static unsigned find_rebin(const TH1 *const)
Definition: LA_Fitter.cc:104
LA_Filler_Fitter::PROB1
Definition: LA_Filler_Fitter.h:53
LA_Filler_Fitter.h
min
T min(T a, T b)
Definition: MathUtil.h:58
LA_Filler_Fitter::rms_profile
static TH1 * rms_profile(const std::string, const TProfile *const)
Definition: LA_Fitter.cc:121
Book::book
TH1 * book(string_t name, TH1 *const hist)
Definition: Book.h:42
hltPixelTracks_cff.chi2
chi2
Definition: hltPixelTracks_cff.py:25
SiStripPI::rms
Definition: SiStripPayloadInspectorHelper.h:169
python.cmstools.all
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
Book
Definition: Book.h:16
LA_Filler_Fitter::AVGV2
Definition: LA_Filler_Fitter.h:54
EDMException.h
alignCSCRings.s
s
Definition: alignCSCRings.py:92
LA_Filler_Fitter::RMSV3
Definition: LA_Filler_Fitter.h:57
SymmetryFit::symmetryChi2
static TH1 * symmetryChi2(std::string, const std::vector< TH1 * > &, const std::pair< unsigned, unsigned >)
Definition: SymmetryFit.cc:6
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
mps_fire.end
end
Definition: mps_fire.py:242
LA_Filler_Fitter::method
static std::string method(Method m, bool fit=true)
Definition: LA_Filler_Fitter.h:61
Book::end
iterator end(string_t re=".*")
Definition: Book.h:65
LA_Filler_Fitter::WIDTH
Definition: LA_Filler_Fitter.h:51
Book::find
iterator find(string_t name, string_t re=".*")
Definition: Book.h:73
Book::erase
void erase(string_t name)
Definition: Book.h:84
LA_Filler_Fitter::make_and_fit_symmchi2
static void make_and_fit_symmchi2(Book &)
Definition: LA_Fitter.cc:40
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
gpuVertexFinder::hist
__shared__ Hist hist
Definition: gpuClusterTracksDBSCAN.h:48
LA_Filler_Fitter::subset_probability
static TH1 * subset_probability(const std::string name, const TH1 *const, const TH1 *const)
Definition: LA_Fitter.cc:135
LA_Filler_Fitter::AVGV3
Definition: LA_Filler_Fitter.h:55
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Book::iterator
boost::filter_iterator< match_name, book_t::iterator > iterator
Definition: Book.h:55
Book::begin
iterator begin(string_t re=".*")
Definition: Book.h:57
TtFullHadDaughter::B
static const std::string B
Definition: TtFullHadronicEvent.h:9
T
long double T
Definition: Basic3DVectorLD.h:48
LA_Filler_Fitter::fit_width_profile
static void fit_width_profile(Book &)
Definition: LA_Fitter.cc:8
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
LA_Filler_Fitter::RMSV2
Definition: LA_Filler_Fitter.h:56
dqmMemoryStats.total
total
Definition: dqmMemoryStats.py:152
trigObjTnPSource_cfi.bins
bins
Definition: trigObjTnPSource_cfi.py:20
mps_check.guess
string guess
Definition: mps_check.py:394
newFWLiteAna.base
base
Definition: newFWLiteAna.py:92
TtFullHadEvtBuilder_cfi.prob
prob
Definition: TtFullHadEvtBuilder_cfi.py:33
fit
Definition: CombinedChiSquaredLikelihood.h:6