35 Double_t Reconstruct( Double_t xt, TRandom3& R )
38 const Double_t cutdummy = -99999.0;
39 Double_t xeff = 0.3 + (1.0 - 0.3)/20.0*(xt + 10.0);
40 Double_t x = R.Rndm();
41 if (x > xeff)
return cutdummy;
43 Double_t xsmear= R.Gaus(-2.5,0.2);
48 void TSVDUnfoldExample()
50 gROOT->SetStyle(
"Plain");
51 gStyle->SetOptStat(0);
55 const Double_t cutdummy= -99999.0;
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++) {
75 Double_t xt = R.BreitWigner(0.3, 2.5);
77 Double_t x = Reconstruct( xt, R );
85 for (Int_t i=0; i<10000; i++) {
86 Double_t xt = R.Gaus(0.0, 2.0);
88 Double_t x = Reconstruct( xt, R );
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++) {
100 statcov->SetBinContent(i,i,data->GetBinError(i)*data->GetBinError(i));
107 TSVDUnfold *tsvdunf =
new TSVDUnfold( data, statcov, bini, xini, Adet );
110 tsvdunf->SetNormalize( kFALSE );
115 TH1D* unfres = tsvdunf->Unfold( 13 );
119 TH1D* ddist = tsvdunf->GetD();
122 TH1D* svdist = tsvdunf->GetSV();
128 TH2D* ustatcov = tsvdunf->GetUnfoldCovMatrix( statcov, 100 );
132 TH2D* uadetcov = tsvdunf->GetAdetCovMatrix( 100 );
135 ustatcov->Add( uadetcov );
138 TH2D* utaucov = tsvdunf->GetXtau();
139 utaucov->Add( uadetcov );
142 TH2D* uinvcov = tsvdunf->GetXinv();
148 for (
int i=1; i<=unfres->GetNbinsX(); i++) {
149 unfres->SetBinError(i, TMath::Sqrt(utaucov->GetBinContent(i,i)));
153 xini->Scale(0.7*datatrue->Integral()/xini->Integral());
155 TLegend *leg =
new TLegend(0.58,0.60,0.99,0.88);
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 );
167 TVirtualPad * c11 = c1->cd(1);
169 TH1D* frame =
new TH1D( *unfres );
170 frame->SetTitle(
"Unfolding toy example with TSVDUnfold" );
171 frame->GetXaxis()->SetTitle(
"x variable" );
172 frame->GetYaxis()->SetTitle(
"Events" );
173 frame->GetXaxis()->SetTitleOffset( 1.25 );
174 frame->GetYaxis()->SetTitleOffset( 1.29 );
177 data->SetLineStyle(2);
178 data->SetLineColor(4);
179 data->SetLineWidth(2);
180 unfres->SetMarkerStyle(20);
181 datatrue->SetLineColor(2);
182 datatrue->SetLineWidth(2);
183 xini->SetLineStyle(2);
184 xini->SetLineColor(8);
185 xini->SetLineWidth(2);
189 unfres->Draw(
"same");
190 datatrue->Draw(
"same");
197 TVirtualPad * c12 = c1->cd(2);
199 TVirtualPad * c2 = c12->cd(1);
200 c2->SetRightMargin ( 0.15 );
202 TH2D* covframe =
new TH2D( *ustatcov );
203 covframe->SetTitle(
"TSVDUnfold covariance matrix" );
204 covframe->GetXaxis()->SetTitle(
"x variable" );
205 covframe->GetYaxis()->SetTitle(
"x variable" );
206 covframe->GetXaxis()->SetTitleOffset( 1.25 );
207 covframe->GetYaxis()->SetTitleOffset( 1.29 );
210 ustatcov->SetLineWidth( 2 );
211 ustatcov->Draw(
"colzsame" );
214 TVirtualPad * c3 = c12->cd(2);
217 TLine *line =
new TLine( 0.,1.,40.,1. );
218 line->SetLineStyle(2);
220 TH1D* dframe =
new TH1D( *ddist );
221 dframe->SetTitle(
"TSVDUnfold |d_{i}|" );
222 dframe->GetXaxis()->SetTitle(
"i" );
223 dframe->GetYaxis()->SetTitle(
"|d_{i}|" );
224 dframe->GetXaxis()->SetTitleOffset( 1.25 );
225 dframe->GetYaxis()->SetTitleOffset( 1.29 );
226 dframe->SetMinimum( 0.001 );
229 ddist->SetLineWidth( 2 );
230 ddist->Draw(
"same" );