39 Double_t xeff = 0.3 + (1.0 - 0.3)/20.0*(xt + 10.0);
41 if (
x > xeff)
return cutdummy;
48void TSVDUnfoldExample()
50 gROOT->SetStyle(
"Plain");
62 TH1D *xini =
new TH1D(
"xini",
"MC truth", nbins, -10.0, 10.0);
63 TH1D *bini =
new TH1D(
"bini",
"MC reco", nbins, -10.0, 10.0);
64 TH2D *Adet =
new TH2D(
"Adet",
"detector response", nbins, -10.0, 10.0, nbins, -10.0, 10.0);
67 TH1D *
data =
new TH1D(
"data",
"data", nbins, -10.0, 10.0);
69 TH1D *datatrue =
new TH1D(
"datatrue",
"data truth", nbins, -10.0, 10.0);
71 TH2D *statcov =
new TH2D(
"statcov",
"covariance matrix", nbins, -10.0, 10.0, nbins, -10.0, 10.0);
74 for (
Int_t i= 0; i<100000; i++) {
85 for (
Int_t i=0; i<10000; i++) {
93 cout <<
"Created toy distributions and errors for: " << endl;
94 cout <<
"... \"true MC\" and \"reconstructed (smeared) MC\"" << endl;
95 cout <<
"... \"true data\" and \"reconstructed (smeared) data\"" << endl;
96 cout <<
"... the \"detector response matrix\"" << endl;
99 for (
int i=1; i<=
data->GetNbinsX(); i++) {
135 ustatcov->
Add( uadetcov );
139 utaucov->
Add( uadetcov );
148 for (
int i=1; i<=unfres->
GetNbinsX(); i++) {
156 leg->SetBorderSize(0);
157 leg->SetFillColor(0);
158 leg->SetFillStyle(0);
159 leg->AddEntry(unfres,
"Unfolded Data",
"p");
160 leg->AddEntry(datatrue,
"True Data",
"l");
161 leg->AddEntry(
data,
"Reconstructed Data",
"l");
162 leg->AddEntry(xini,
"True MC",
"l");
164 TCanvas *
c1 =
new TCanvas(
"c1",
"Unfolding toy example with TSVDUnfold", 1000, 900 );
170 frame->
SetTitle(
"Unfolding toy example with TSVDUnfold" );
177 data->SetLineStyle(2);
178 data->SetLineColor(4);
179 data->SetLineWidth(2);
189 unfres->
Draw(
"same");
190 datatrue->
Draw(
"same");
202 TH2D* covframe =
new TH2D( *ustatcov );
203 covframe->
SetTitle(
"TSVDUnfold covariance matrix" );
211 ustatcov->
Draw(
"colzsame" );
221 dframe->
SetTitle(
"TSVDUnfold |d_{i}|" );
230 ddist->
Draw(
"same" );
#define R(a, b, c, d, e, f, g, h, i)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
R__EXTERN TStyle * gStyle
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
virtual void SetLineColor(Color_t lcolor)
Set the line color.
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
virtual void SetRightMargin(Float_t rightmargin)
Set Pad right margin in fraction of the pad width.
1-D histogram with a double per channel (see TH1 documentation)
void SetTitle(const char *title) override
Change/set the title.
virtual Int_t GetNbinsX() const
virtual Bool_t Add(TF1 *h1, Double_t c1=1, Option_t *option="")
Performs the operation: this = this + c1*f1 if errors are defined (see TH1::Sumw2),...
virtual void SetBinError(Int_t bin, Double_t error)
Set the bin Error Note that this resets the bin eror option to be of Normal Type and for the non-empt...
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
void Draw(Option_t *option="") override
Draw this histogram with options.
virtual void SetMinimum(Double_t minimum=-1111)
virtual Double_t Integral(Option_t *option="") const
Return integral of bin contents.
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
2-D histogram with a double per channel (see TH1 documentation)
void SetBinContent(Int_t bin, Double_t content) override
Set bin content.
Double_t GetBinContent(Int_t binx, Int_t biny) const override
Int_t Fill(Double_t) override
Invalid Fill method.
This class displays a legend box (TPaveText) containing several legend entries.
Use the TLine constructor to create a simple line.
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Random number generator class based on M.
SVD Approach to Data Unfolding.
TH1D * GetSV() const
Returns singular values vector.
TH2D * GetXtau() const
Returns the computed regularized covariance matrix corresponding to total uncertainties on measured s...
TH1D * Unfold(Int_t kreg)
Perform the unfolding with regularisation parameter kreg.
void SetNormalize(Bool_t normalize)
TH1D * GetD() const
Returns d vector (for choosing appropriate regularisation)
TH2D * GetUnfoldCovMatrix(const TH2D *cov, Int_t ntoys, Int_t seed=1)
Determine for given input error matrix covariance matrix of unfolded spectrum from toy simulation giv...
TH2D * GetAdetCovMatrix(Int_t ntoys, Int_t seed=1)
Determine covariance matrix of unfolded spectrum from finite statistics in response matrix using pseu...
TH2D * GetXinv() const
Returns the computed inverse of the covariance matrix.
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
TVirtualPad is an abstract base class for the Pad and Canvas classes.
virtual TVirtualPad * cd(Int_t subpadnumber=0)=0
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)=0
virtual void SetLogy(Int_t value=1)=0
Double_t Sqrt(Double_t x)
Returns the square root of x.