Logo ROOT  
Reference Guide
testUnfold1.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_unfold
3/// \notebook
4/// Test program for the classes TUnfold and related.
5///
6/// 1. Generate Monte Carlo and Data events
7/// The events consist of
8/// signal
9/// background
10///
11/// The signal is a resonance. It is generated with a Breit-Wigner,
12/// smeared by a Gaussian
13///
14/// 2. Unfold the data. The result is:
15/// The background level
16/// The shape of the resonance, corrected for detector effects
17///
18/// Systematic errors from the MC shape variation are included
19/// and propagated to the result
20///
21/// 3. fit the unfolded distribution, including the correlation matrix
22///
23/// 4. save six plots to a file `testUnfold1.ps`
24/// ~~~
25/// 1 2 3
26/// 4 5 6
27/// ~~~
28/// 1. 2d-plot of the matrix describing the migrations
29/// 2. generator-level distributions
30/// - blue: unfolded data, total errors
31/// - green: unfolded data, statistical errors
32/// - red: generated data
33/// - black: fit to green data points
34/// 3. detector level distributions
35/// - blue: unfolded data, folded back through the matrix
36/// - black: Monte Carlo (with wrong peal position)
37/// - blue: data
38/// 4. global correlation coefficients
39/// 5. \f$ \chi^2 \f$ as a function of \f$ log(\tau) \f$
40/// the star indicates the final choice of \f$ \tau \f$
41/// 6. the L curve
42///
43/// \macro_output
44/// \macro_image
45/// \macro_code
46///
47/// **Version 17.6, in parallel to changes in TUnfold**
48///
49/// #### History:
50/// - Version 17.5, in parallel to changes in TUnfold
51/// - Version 17.4, in parallel to changes in TUnfold
52/// - Version 17.3, in parallel to changes in TUnfold
53/// - Version 17.2, in parallel to changes in TUnfold
54/// - Version 17.1, in parallel to changes in TUnfold
55/// - Version 17.0, updated for using the classes TUnfoldDensity, TUnfoldBinning
56/// - Version 16.1, parallel to changes in TUnfold
57/// - Version 16.0, parallel to changes in TUnfold
58/// - Version 15, with automated L-curve scan
59/// - Version 14, with changes in TUnfoldSys.cxx
60/// - Version 13, include test of systematic errors
61/// - Version 12, catch error when defining the input
62/// - Version 11, print chi**2 and number of degrees of freedom
63/// - Version 10, with bug-fix in TUnfold.cxx
64/// - Version 9, with bug-fix in TUnfold.cxx and TUnfold.h
65/// - Version 8, with bug-fix in TUnfold.cxx and TUnfold.h
66/// - Version 7, with bug-fix in TUnfold.cxx and TUnfold.h
67/// - Version 6a, fix problem with dynamic array allocation under windows
68/// - Version 6, bug-fixes in TUnfold.C
69/// - Version 5, replace main() by testUnfold1()
70/// - Version 4, with bug-fix in TUnfold.C
71/// - Version 3, with bug-fix in TUnfold.C
72/// - Version 2, with changed ScanLcurve() arguments
73/// - Version 1, remove L curve analysis, use ScanLcurve() method instead
74/// - Version 0, L curve analysis included here
75///
76/// This file is part of TUnfold.
77///
78/// TUnfold is free software: you can redistribute it and/or modify
79/// it under the terms of the GNU General Public License as published by
80/// the Free Software Foundation, either version 3 of the License, or
81/// (at your option) any later version.
82///
83/// TUnfold is distributed in the hope that it will be useful,
84/// but WITHOUT ANY WARRANTY; without even the implied warranty of
85/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
86/// GNU General Public License for more details.
87///
88/// You should have received a copy of the GNU General Public License
89/// along with TUnfold. If not, see <http://www.gnu.org/licenses/>.
90///
91/// \author Stefan Schmitt DESY, 14.10.2008
92
93
94#include <TError.h>
95#include <TMath.h>
96#include <TCanvas.h>
97#include <TRandom3.h>
98#include <TFitter.h>
99#include <TF1.h>
100#include <TStyle.h>
101#include <TVector.h>
102#include <TGraph.h>
103
104#include "TUnfoldDensity.h"
105
106// #define VERBOSE_LCURVE_SCAN
107
108using namespace std;
109
110TRandom *rnd=0;
111
112TH2 *gHistInvEMatrix;
113
114TVirtualFitter *gFitter=0;
115
116void chisquare_corr(Int_t &npar, Double_t * /*gin */, Double_t &f, Double_t *u, Int_t /*flag */) {
117 // Minimization function for H1s using a Chisquare method
118 // only one-dimensional histograms are supported
119 // Correlated errors are taken from an external inverse covariance matrix
120 // stored in a 2-dimensional histogram
121 Double_t x;
122 TH1 *hfit = (TH1*)gFitter->GetObjectFit();
123 TF1 *f1 = (TF1*)gFitter->GetUserFunc();
124
125 f1->InitArgs(&x,u);
126 npar = f1->GetNpar();
127 f = 0;
128
129 Int_t npfit = 0;
130 Int_t nPoints=hfit->GetNbinsX();
131 Double_t *df=new Double_t[nPoints];
132 for (Int_t i=0;i<nPoints;i++) {
133 x = hfit->GetBinCenter(i+1);
135 df[i] = f1->EvalPar(&x,u)-hfit->GetBinContent(i+1);
136 if (TF1::RejectedPoint()) df[i]=0.0;
137 else npfit++;
138 }
139 for (Int_t i=0;i<nPoints;i++) {
140 for (Int_t j=0;j<nPoints;j++) {
141 f += df[i]*df[j]*gHistInvEMatrix->GetBinContent(i+1,j+1);
142 }
143 }
144 delete[] df;
145 f1->SetNumberFitPoints(npfit);
146}
147
148Double_t bw_func(Double_t *x,Double_t *par) {
149 Double_t dm=x[0]-par[1];
150 return par[0]/(dm*dm+par[2]*par[2]);
151}
152
153
154// generate an event
155// output:
156// negative mass: background event
157// positive mass: signal event
158Double_t GenerateEvent(Double_t bgr, // relative fraction of background
159 Double_t mass, // peak position
160 Double_t gamma) // peak width
161{
162 Double_t t;
163 if(rnd->Rndm()>bgr) {
164 // generate signal event
165 // with positive mass
166 do {
167 do {
168 t=rnd->Rndm();
169 } while(t>=1.0);
170 t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass;
171 } while(t<=0.0);
172 return t;
173 } else {
174 // generate background event
175 // generate events following a power-law distribution
176 // f(E) = K * TMath::power((E0+E),N0)
177 static Double_t const E0=2.4;
178 static Double_t const N0=2.9;
179 do {
180 do {
181 t=rnd->Rndm();
182 } while(t>=1.0);
183 // the mass is returned negative
184 // In our example a convenient way to indicate it is a background event.
185 t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0;
186 } while(t>=0.0);
187 return t;
188 }
189}
190
191// smear the event to detector level
192// input:
193// mass on generator level (mTrue>0 !)
194// output:
195// mass on detector level
196Double_t DetectorEvent(Double_t mTrue) {
197 // smear by double-gaussian
198 static Double_t frac=0.1;
199 static Double_t wideBias=0.03;
200 static Double_t wideSigma=0.5;
201 static Double_t smallBias=0.0;
202 static Double_t smallSigma=0.1;
203 if(rnd->Rndm()>frac) {
204 return rnd->Gaus(mTrue+smallBias,smallSigma);
205 } else {
206 return rnd->Gaus(mTrue+wideBias,wideSigma);
207 }
208}
209
210int testUnfold1()
211{
212 // switch on histogram errors
214
215 // show fit result
216 gStyle->SetOptFit(1111);
217
218 // random generator
219 rnd=new TRandom3();
220
221 // data and MC luminosity, cross-section
222 Double_t const luminosityData=100000;
223 Double_t const luminosityMC=1000000;
224 Double_t const crossSection=1.0;
225
226 Int_t const nDet=250;
227 Int_t const nGen=100;
228 Double_t const xminDet=0.0;
229 Double_t const xmaxDet=10.0;
230 Double_t const xminGen=0.0;
231 Double_t const xmaxGen=10.0;
232
233 //============================================
234 // generate MC distribution
235 //
236 TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen);
237 TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet);
238 TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",
239 nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
240 Int_t neventMC=rnd->Poisson(luminosityMC*crossSection);
241 for(Int_t i=0;i<neventMC;i++) {
242 Double_t mGen=GenerateEvent(0.3, // relative fraction of background
243 4.0, // peak position in MC
244 0.2); // peak width in MC
245 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
246 // the generated mass is negative for background
247 // and positive for signal
248 // so it will be filled in the underflow bin
249 // this is very convenient for the unfolding:
250 // the unfolded result will contain the number of background
251 // events in the underflow bin
252
253 // generated MC distribution (for comparison only)
254 histMgenMC->Fill(mGen,luminosityData/luminosityMC);
255 // reconstructed MC distribution (for comparison only)
256 histMdetMC->Fill(mDet,luminosityData/luminosityMC);
257
258 // matrix describing how the generator input migrates to the
259 // reconstructed level. Unfolding input.
260 // NOTE on underflow/overflow bins:
261 // (1) the detector level under/overflow bins are used for
262 // normalisation ("efficiency" correction)
263 // in our toy example, these bins are populated from tails
264 // of the initial MC distribution.
265 // (2) the generator level underflow/overflow bins are
266 // unfolded. In this example:
267 // underflow bin: background events reconstructed in the detector
268 // overflow bin: signal events generated at masses > xmaxDet
269 // for the unfolded result these bins will be filled
270 // -> the background normalisation will be contained in the underflow bin
271 histMdetGenMC->Fill(mDet,mGen,luminosityData/luminosityMC);
272 }
273
274 //============================================
275 // generate alternative MC
276 // this will be used to derive a systematic error due to MC
277 // parameter uncertainties
278 TH2D *histMdetGenSysMC=new TH2D("MdetgenSysMC",";mass(det);mass(gen)",
279 nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
280 neventMC=rnd->Poisson(luminosityMC*crossSection);
281 for(Int_t i=0;i<neventMC;i++) {
282 Double_t mGen=GenerateEvent
283 (0.5, // relative fraction of background
284 3.6, // peak position in MC with systematic shift
285 0.15); // peak width in MC
286 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
287 histMdetGenSysMC->Fill(mDet,mGen,luminosityData/luminosityMC);
288 }
289
290 //============================================
291 // generate data distribution
292 //
293 TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen);
294 TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet);
295 Int_t neventData=rnd->Poisson(luminosityData*crossSection);
296 for(Int_t i=0;i<neventData;i++) {
297 Double_t mGen=GenerateEvent(0.4, // relative fraction of background
298 3.8, // peak position in data
299 0.15); // peak width in data
300 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
301 // generated data mass for comparison plots
302 // for real data, we do not have this histogram
303 histMgenData->Fill(mGen);
304
305 // reconstructed mass, unfolding input
306 histMdetData->Fill(mDet);
307 }
308
309 //=========================================================================
310 // divide by bin width to get density distributions
311 TH1D *histDensityGenData=new TH1D("DensityGenData",";mass(gen)",
312 nGen,xminGen,xmaxGen);
313 TH1D *histDensityGenMC=new TH1D("DensityGenMC",";mass(gen)",
314 nGen,xminGen,xmaxGen);
315 for(Int_t i=1;i<=nGen;i++) {
316 histDensityGenData->SetBinContent(i,histMgenData->GetBinContent(i)/
317 histMgenData->GetBinWidth(i));
318 histDensityGenMC->SetBinContent(i,histMgenMC->GetBinContent(i)/
319 histMgenMC->GetBinWidth(i));
320 }
321
322 //=========================================================================
323 // set up the unfolding
324 // define migration matrix
325 TUnfoldDensity unfold(histMdetGenMC,TUnfold::kHistMapOutputVert);
326
327 // define input and bias scheme
328 // do not use the bias, because MC peak may be at the wrong place
329 // watch out for error codes returned by the SetInput method
330 // errors larger or equal 10000 are fatal:
331 // the data points specified as input are not sufficient to constrain the
332 // unfolding process
333 if(unfold.SetInput(histMdetData)>=10000) {
334 std::cout<<"Unfolding result may be wrong\n";
335 }
336
337 //========================================================================
338 // the unfolding is done here
339 //
340 // scan L curve and find best point
341 Int_t nScan=30;
342 // use automatic L-curve scan: start with taumin=taumax=0.0
343 Double_t tauMin=0.0;
344 Double_t tauMax=0.0;
345 Int_t iBest;
346 TSpline *logTauX,*logTauY;
347 TGraph *lCurve;
348
349 // if required, report Info messages (for debugging the L-curve scan)
350#ifdef VERBOSE_LCURVE_SCAN
351 Int_t oldinfo=gErrorIgnoreLevel;
353#endif
354 // this method scans the parameter tau and finds the kink in the L curve
355 // finally, the unfolding is done for the best choice of tau
356 iBest=unfold.ScanLcurve(nScan,tauMin,tauMax,&lCurve,&logTauX,&logTauY);
357
358 // if required, switch to previous log-level
359#ifdef VERBOSE_LCURVE_SCAN
360 gErrorIgnoreLevel=oldinfo;
361#endif
362
363 //==========================================================================
364 // define a correlated systematic error
365 // for example, assume there is a 10% correlated error for all reconstructed
366 // masses larger than 7
367 Double_t SYS_ERROR1_MSTART=6;
368 Double_t SYS_ERROR1_SIZE=0.1;
369 TH2D *histMdetGenSys1=new TH2D("Mdetgensys1",";mass(det);mass(gen)",
370 nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
371 for(Int_t i=0;i<=nDet+1;i++) {
372 if(histMdetData->GetBinCenter(i)>=SYS_ERROR1_MSTART) {
373 for(Int_t j=0;j<=nGen+1;j++) {
374 histMdetGenSys1->SetBinContent(i,j,SYS_ERROR1_SIZE);
375 }
376 }
377 }
378 unfold.AddSysError(histMdetGenSysMC,"SYSERROR_MC",TUnfold::kHistMapOutputVert,
380 unfold.AddSysError(histMdetGenSys1,"SYSERROR1",TUnfold::kHistMapOutputVert,
382
383 //==========================================================================
384 // print some results
385 //
386 std::cout<<"tau="<<unfold.GetTau()<<"\n";
387 std::cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
388 <<" / "<<unfold.GetNdf()<<"\n";
389 std::cout<<"chi**2(sys)="<<unfold.GetChi2Sys()<<"\n";
390
391
392 //==========================================================================
393 // create graphs with one point to visualize the best choice of tau
394 //
395 Double_t t[1],x[1],y[1];
396 logTauX->GetKnot(iBest,t[0],x[0]);
397 logTauY->GetKnot(iBest,t[0],y[0]);
398 TGraph *bestLcurve=new TGraph(1,x,y);
399 TGraph *bestLogTauLogChi2=new TGraph(1,t,x);
400
401 //==========================================================================
402 // retrieve results into histograms
403
404 // get unfolded distribution
405 TH1 *histMunfold=unfold.GetOutput("Unfolded");
406
407 // get unfolding result, folded back
408 TH1 *histMdetFold=unfold.GetFoldedOutput("FoldedBack");
409
410 // get error matrix (input distribution [stat] errors only)
411 // TH2D *histEmatData=unfold.GetEmatrix("EmatData");
412
413 // get total error matrix:
414 // migration matrix uncorrelated and correlated systematic errors
415 // added in quadrature to the data statistical errors
416 TH2 *histEmatTotal=unfold.GetEmatrixTotal("EmatTotal");
417
418 // create data histogram with the total errors
419 TH1D *histTotalError=
420 new TH1D("TotalError",";mass(gen)",nGen,xminGen,xmaxGen);
421 for(Int_t bin=1;bin<=nGen;bin++) {
422 histTotalError->SetBinContent(bin,histMunfold->GetBinContent(bin));
423 histTotalError->SetBinError
424 (bin,TMath::Sqrt(histEmatTotal->GetBinContent(bin,bin)));
425 }
426
427 // get global correlation coefficients
428 // for this calculation one has to specify whether the
429 // underflow/overflow bins are included or not
430 // default: include all bins
431 // here: exclude underflow and overflow bins
432 TH1 *histRhoi=unfold.GetRhoItotal("rho_I",
433 0, // use default title
434 0, // all distributions
435 "*[UO]", // discard underflow and overflow bins on all axes
436 kTRUE, // use original binning
437 &gHistInvEMatrix // store inverse of error matrix
438 );
439
440 //======================================================================
441 // fit Breit-Wigner shape to unfolded data, using the full error matrix
442 // here we use a "user" chi**2 function to take into account
443 // the full covariance matrix
444
445 gFitter=TVirtualFitter::Fitter(histMunfold);
446 gFitter->SetFCN(chisquare_corr);
447
448 TF1 *bw=new TF1("bw",bw_func,xminGen,xmaxGen,3);
449 bw->SetParameter(0,1000.);
450 bw->SetParameter(1,3.8);
451 bw->SetParameter(2,0.2);
452
453 // for (wrong!) fitting without correlations, drop the option "U"
454 // here.
455 histMunfold->Fit(bw,"UE");
456
457 //=====================================================================
458 // plot some histograms
460 output.Divide(3,2);
461
462 // Show the matrix which connects input and output
463 // There are overflow bins at the bottom, not shown in the plot
464 // These contain the background shape.
465 // The overflow bins to the left and right contain
466 // events which are not reconstructed. These are necessary for proper MC
467 // normalisation
468 output.cd(1);
469 histMdetGenMC->Draw("BOX");
470
471 // draw generator-level distribution:
472 // data (red) [for real data this is not available]
473 // MC input (black) [with completely wrong peak position and shape]
474 // unfolded data (blue)
475 output.cd(2);
476 histTotalError->SetLineColor(kBlue);
477 histTotalError->Draw("E");
478 histMunfold->SetLineColor(kGreen);
479 histMunfold->Draw("SAME E1");
480 histDensityGenData->SetLineColor(kRed);
481 histDensityGenData->Draw("SAME");
482 histDensityGenMC->Draw("SAME HIST");
483
484 // show detector level distributions
485 // data (red)
486 // MC (black) [with completely wrong peak position and shape]
487 // unfolded data (blue)
488 output.cd(3);
489 histMdetFold->SetLineColor(kBlue);
490 histMdetFold->Draw();
491 histMdetMC->Draw("SAME HIST");
492
493 TH1 *histInput=unfold.GetInput("Minput",";mass(det)");
494
495 histInput->SetLineColor(kRed);
496 histInput->Draw("SAME");
497
498 // show correlation coefficients
499 output.cd(4);
500 histRhoi->Draw();
501
502 // show tau as a function of chi**2
503 output.cd(5);
504 logTauX->Draw();
505 bestLogTauLogChi2->SetMarkerColor(kRed);
506 bestLogTauLogChi2->Draw("*");
507
508 // show the L curve
509 output.cd(6);
510 lCurve->Draw("AL");
511 bestLcurve->SetMarkerColor(kRed);
512 bestLcurve->Draw("*");
513
514 output.SaveAs("testUnfold1.ps");
515
516 return 0;
517}
#define f(i)
Definition: RSha256.hxx:104
const Bool_t kFALSE
Definition: RtypesCore.h:90
const Bool_t kTRUE
Definition: RtypesCore.h:89
@ kRed
Definition: Rtypes.h:64
@ kGreen
Definition: Rtypes.h:64
@ kBlue
Definition: Rtypes.h:64
R__EXTERN Int_t gErrorIgnoreLevel
Definition: TError.h:105
const Int_t kInfo
Definition: TError.h:37
R__EXTERN TStyle * gStyle
Definition: TStyle.h:410
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
The Canvas class.
Definition: TCanvas.h:27
1-Dim function class
Definition: TF1.h:210
static void RejectPoint(Bool_t reject=kTRUE)
Static function to set the global flag to reject points the fgRejectPoint global flag is tested by al...
Definition: TF1.cxx:3647
virtual Int_t GetNpar() const
Definition: TF1.h:475
virtual void SetNumberFitPoints(Int_t npfits)
Definition: TF1.h:618
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=0)
Evaluate function with given coordinates and parameters.
Definition: TF1.cxx:1461
virtual void InitArgs(const Double_t *x, const Double_t *params)
Initialize parameters addresses.
Definition: TF1.cxx:2457
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition: TF1.cxx:3656
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:628
A TGraph is an object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:760
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:614
The TH1 histogram class.
Definition: TH1.h:56
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
Definition: TH1.cxx:8597
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3808
virtual Int_t GetNbinsX() const
Definition: TH1.h:292
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...
Definition: TH1.cxx:8662
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3275
static void SetDefaultSumw2(Bool_t sumw2=kTRUE)
When this static function is called with sumw2=kTRUE, all new histograms will automatically activate ...
Definition: TH1.cxx:6330
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition: TH1.cxx:8678
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4907
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width for 1D histogram.
Definition: TH1.cxx:8619
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:292
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:294
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:88
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2480
Random number generator class based on M.
Definition: TRandom3.h:27
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:27
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:263
virtual Int_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
Definition: TRandom.cxx:391
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom.cxx:541
Base class for spline implementation containing the Draw/Paint methods.
Definition: TSpline.h:22
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TSpline.cxx:97
virtual void GetKnot(Int_t i, Double_t &x, Double_t &y) const =0
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition: TStyle.cxx:1542
An algorithm to unfold distributions from detector to truth level.
@ kSysErrModeRelative
matrix gives the relative shifts
Definition: TUnfoldSys.h:108
@ kSysErrModeMatrix
matrix is an alternative to the default matrix, the errors are the difference to the original matrix
Definition: TUnfoldSys.h:104
@ kHistMapOutputVert
truth level on y-axis of the response matrix
Definition: TUnfold.h:145
Abstract Base Class for Fitting.
virtual TObject * GetObjectFit() const
virtual void SetFCN(void(*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t))
To set the address of the minimization objective function called by the native compiler (see function...
virtual TObject * GetUserFunc() const
static TVirtualFitter * Fitter(TObject *obj, Int_t maxpar=25)
Static function returning a pointer to the current fitter.
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
TF1 * f1
Definition: legend1.C:11
double gamma(double x)
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:725
constexpr Double_t Pi()
Definition: TMath.h:38
Double_t Tan(Double_t)
Definition: TMath.h:635
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
static void output(int code)
Definition: gifencode.c:226