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 
10  for(Book::iterator it = book.begin(".*"+method(WIDTH)); it!=book.end(); ++it) {
11  it->second->SetTitle("Mean Cluster Width;tan#theta_{t}");
12  TH1* const p = it->second;
13  if(p->GetEntries() < 400) { delete p; book[it->first]=nullptr; continue;}
14  p->SetTitle(";tan#theta_{t};");
15  const float min = p->GetMinimum();
16  const float max = p->GetMaximum();
17  float xofmin = p->GetBinCenter(p->GetMinimumBin()); if( xofmin>0.0 || xofmin<-0.15) xofmin = -0.05;
18  const float xofmax = p->GetBinCenter(p->GetMaximumBin());
19 
20  TF1* const fit = new TF1("LA_profile_fit","[2]*(TMath::Abs(x-[0]))+[1]",-1,1);
21  fit->SetParLimits(0,-0.15,0.01);
22  fit->SetParLimits(1, 0.6*min, 1.25*min );
23  fit->SetParLimits(2,0.1,10);
24  fit->SetParameters( xofmin, min, (max-min) / fabs( xofmax - xofmin ) );
25 
26  int badfit = p->Fit(fit,"IEQ","",-.5,.3);
27  if( badfit ) badfit = p->Fit(fit,"IEQ","", -.46,.26);
28  if( badfit ) { book.erase(it); }
29  }
30 }
31 
34  for(Book::iterator it = book.begin(".*_all"); it!=book.end(); ++it) {
35  const std::string base = boost::erase_all_copy(it->first,"_all");
36 
37  std::vector<Book::iterator> rebin_hists;
38  const Book::iterator& all = it; rebin_hists.push_back(all);
39  Book::iterator w1 = book.find(base+"_w1"); rebin_hists.push_back(w1);
40  Book::iterator var_w2 = book.find(base+method(AVGV2,false)); rebin_hists.push_back(var_w2);
41  Book::iterator var_w3 = book.find(base+method(AVGV3,false)); rebin_hists.push_back(var_w3);
42 
43  const unsigned rebin = std::max( var_w2==book.end() ? 0 : find_rebin(var_w2->second),
44  var_w3==book.end() ? 0 : find_rebin(var_w3->second) );
45  for(const auto& it : rebin_hists) if(it!=book.end()) it->second->Rebin( rebin>1 ? rebin<7 ? rebin : 6 : 1);
46 
47  TH1* const prob_w1 = w1==book.end() ? nullptr : subset_probability( base+method(PROB1,false) ,w1->second,all->second);
48  TH1* const rmsv_w2 = var_w2==book.end() ? nullptr : rms_profile( base+method(RMSV2,false), (TProfile*const)var_w2->second);
49  TH1* const rmsv_w3 = var_w3==book.end() ? nullptr : rms_profile( base+method(RMSV3,false), (TProfile*const)var_w3->second);
50 
51  std::vector<TH1*> fit_hists;
52  if(prob_w1) {
53  book.book(base+method(PROB1,false),prob_w1);
54  fit_hists.push_back(prob_w1); prob_w1->SetTitle("Width==1 Probability;tan#theta_{t}-(dx/dz)_{reco}");
55  }
56  if(var_w2!=book.end()) {
57  book.book(base+method(RMSV2,false),rmsv_w2);
58  fit_hists.push_back(var_w2->second); var_w2->second->SetTitle("Width==2 Mean Variance;tan#theta_{t}-(dx/dz)_{reco}");
59  fit_hists.push_back(rmsv_w2); rmsv_w2->SetTitle("Width==2 RMS Variance;tan#theta_{t}-(dx/dz)_{reco}");
60  }
61  if(var_w3!=book.end()) {
62  book.book(base+method(RMSV3,false),rmsv_w3);
63  fit_hists.push_back(var_w3->second); var_w3->second->SetTitle("Width==3 Mean Variance;tan#theta_{t}-(dx/dz)_{reco}");
64  fit_hists.push_back(rmsv_w3); rmsv_w3->SetTitle("Width==3 RMS Variance;tan#theta_{t}-(dx/dz)_{reco}");
65  }
66 
67  if(fit_hists.empty()) continue;
68  const unsigned bins = fit_hists[0]->GetNbinsX();
69  const unsigned guess = fit_hists[0]->FindBin(0);
70  const std::pair<unsigned,unsigned> range(guess-bins/30,guess+bins/30);
71 
72  for(auto const& hist : fit_hists) {
73  TH1*const chi2 = SymmetryFit::symmetryChi2(hist,range);
74  if(chi2) {book.book(chi2->GetName(),chi2); chi2->SetTitle("Symmetry #chi^{2};tan#theta_{t}-(dx/dz)_{reco}");}
75  }
76  }
77 }
78 
79 unsigned LA_Filler_Fitter::
80 find_rebin(const TH1* const hist) {
81  const double mean = hist->GetMean();
82  const double rms = hist->GetRMS();
83  const int begin = std::min( 1, hist->GetXaxis()->FindFixBin(mean-rms));
84  const int end = std::max(hist->GetNbinsX(), hist->GetXaxis()->FindFixBin(mean+rms)) + 1;
85  unsigned current_hole(0), max_hole(0);
86  for(int i=begin; i<end; i++) {
87  if(!hist->GetBinError(i)) current_hole++;
88  else if(current_hole) {max_hole = std::max(current_hole,max_hole); current_hole=0;}
89  }
90  return max_hole+1;
91 }
92 
94 rms_profile(const std::string name, const TProfile* const prof) {
95  const int bins = prof->GetNbinsX();
96  TH1* const rms = new TH1F(name.c_str(),"",bins, prof->GetBinLowEdge(1), prof->GetBinLowEdge(bins) + prof->GetBinWidth(bins) );
97  for(int i = 1; i<=bins; i++) {
98  const double Me = prof->GetBinError(i);
99  const double neff = prof->GetBinEntries(i); //Should be prof->GetBinEffectiveEntries(i);, not availible this version ROOT. This is only ok for unweighted fills
100  rms->SetBinContent(i, Me*sqrt(neff) );
101  rms->SetBinError(i, Me/sqrt(2) );
102  }
103  return rms;
104 }
105 
107 subset_probability(const std::string name, const TH1* const subset, const TH1* const total) {
108  const int bins = subset->GetNbinsX();
109  TH1* const prob = new TH1F(name.c_str(),"",bins, subset->GetBinLowEdge(1), subset->GetBinLowEdge(bins) + subset->GetBinWidth(bins) );
110  for(int i = 1; i<=bins; i++) {
111  const double s = subset->GetBinContent(i);
112  const double T = total->GetBinContent(i);
113  const double B = T-s;
114 
115  const double p = T? s/T : 0;
116  const double perr = T? ( (s&&B)? sqrt(s*s*B+B*B*s)/(T*T) : 1/T ) : 0;
117 
118  prob->SetBinContent(i,p);
119  prob->SetBinError(i,perr);
120  }
121  return prob;
122 }
static TH1 * subset_probability(const std::string name, const TH1 *const , const TH1 *const )
Definition: LA_Fitter.cc:107
static void fit_width_profile(Book &)
Definition: LA_Fitter.cc:9
static TH1 * symmetryChi2(std::string, const std::vector< TH1 * > &, const std::pair< unsigned, unsigned >)
Definition: SymmetryFit.cc:6
T sqrt(T t)
Definition: SSEVec.h:18
iterator begin(string_t re=".*")
Definition: Book.h:48
#define end
Definition: vmac.h:39
T min(T a, T b)
Definition: MathUtil.h:58
base
Make Sure CMSSW is Setup ##.
static const std::string B
static std::string method(Method m, bool fit=true)
static unsigned find_rebin(const TH1 *const)
Definition: LA_Fitter.cc:80
string guess
Definition: mps_check.py:388
static TH1 * rms_profile(const std::string, const TProfile *const)
Definition: LA_Fitter.cc:94
static void make_and_fit_symmchi2(Book &)
Definition: LA_Fitter.cc:33
#define begin
Definition: vmac.h:32
Definition: Book.h:16
long double T
iterator find(string_t name, string_t re=".*")
Definition: Book.h:52
void erase(string_t name)
Definition: Book.h:57
boost::filter_iterator< match_name, book_t::iterator > iterator
Definition: Book.h:46
iterator end(string_t re=".*")
Definition: Book.h:50
TH1 * book(string_t name, TH1 *const hist)
Definition: Book.h:42