Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
108TRandom *rnd=nullptr;
109
110TH2 *gHistInvEMatrix;
111
112TVirtualFitter *gFitter=nullptr;
113
114void chisquare_corr(Int_t &npar, Double_t * /*gin */, Double_t &f, Double_t *u, Int_t /*flag */) {
115 // Minimization function for H1s using a Chisquare method
116 // only one-dimensional histograms are supported
117 // Correlated errors are taken from an external inverse covariance matrix
118 // stored in a 2-dimensional histogram
119 Double_t x;
120 TH1 *hfit = (TH1*)gFitter->GetObjectFit();
121 TF1 *f1 = (TF1*)gFitter->GetUserFunc();
122
123 f1->InitArgs(&x,u);
124 npar = f1->GetNpar();
125 f = 0;
126
127 Int_t npfit = 0;
128 Int_t nPoints=hfit->GetNbinsX();
129 Double_t *df=new Double_t[nPoints];
130 for (Int_t i=0;i<nPoints;i++) {
131 x = hfit->GetBinCenter(i+1);
133 df[i] = f1->EvalPar(&x,u)-hfit->GetBinContent(i+1);
134 if (TF1::RejectedPoint()) df[i]=0.0;
135 else npfit++;
136 }
137 for (Int_t i=0;i<nPoints;i++) {
138 for (Int_t j=0;j<nPoints;j++) {
139 f += df[i]*df[j]*gHistInvEMatrix->GetBinContent(i+1,j+1);
140 }
141 }
142 delete[] df;
143 f1->SetNumberFitPoints(npfit);
144}
145
146Double_t bw_func(Double_t *x,Double_t *par) {
147 Double_t dm=x[0]-par[1];
148 return par[0]/(dm*dm+par[2]*par[2]);
149}
150
151
152// generate an event
153// output:
154// negative mass: background event
155// positive mass: signal event
156Double_t GenerateEvent(Double_t bgr, // relative fraction of background
157 Double_t mass, // peak position
158 Double_t gamma) // peak width
159{
160 Double_t t;
161 if(rnd->Rndm()>bgr) {
162 // generate signal event
163 // with positive mass
164 do {
165 do {
166 t=rnd->Rndm();
167 } while(t>=1.0);
168 t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass;
169 } while(t<=0.0);
170 return t;
171 } else {
172 // generate background event
173 // generate events following a power-law distribution
174 // f(E) = K * TMath::power((E0+E),N0)
175 static Double_t const E0=2.4;
176 static Double_t const N0=2.9;
177 do {
178 do {
179 t=rnd->Rndm();
180 } while(t>=1.0);
181 // the mass is returned negative
182 // In our example a convenient way to indicate it is a background event.
183 t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0;
184 } while(t>=0.0);
185 return t;
186 }
187}
188
189// smear the event to detector level
190// input:
191// mass on generator level (mTrue>0 !)
192// output:
193// mass on detector level
194Double_t DetectorEvent(Double_t mTrue) {
195 // smear by double-gaussian
196 static Double_t frac=0.1;
197 static Double_t wideBias=0.03;
198 static Double_t wideSigma=0.5;
199 static Double_t smallBias=0.0;
200 static Double_t smallSigma=0.1;
201 if(rnd->Rndm()>frac) {
202 return rnd->Gaus(mTrue+smallBias,smallSigma);
203 } else {
204 return rnd->Gaus(mTrue+wideBias,wideSigma);
205 }
206}
207
208int testUnfold1()
209{
210 // switch on histogram errors
212
213 // show fit result
214 gStyle->SetOptFit(1111);
215
216 // random generator
217 rnd=new TRandom3();
218
219 // data and MC luminosity, cross-section
220 Double_t const luminosityData=100000;
221 Double_t const luminosityMC=1000000;
222 Double_t const crossSection=1.0;
223
224 Int_t const nDet=250;
225 Int_t const nGen=100;
226 Double_t const xminDet=0.0;
227 Double_t const xmaxDet=10.0;
228 Double_t const xminGen=0.0;
229 Double_t const xmaxGen=10.0;
230
231 //============================================
232 // generate MC distribution
233 //
234 TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen);
235 TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet);
236 TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",
237 nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
238 Int_t neventMC=rnd->Poisson(luminosityMC*crossSection);
239 for(Int_t i=0;i<neventMC;i++) {
240 Double_t mGen=GenerateEvent(0.3, // relative fraction of background
241 4.0, // peak position in MC
242 0.2); // peak width in MC
243 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
244 // the generated mass is negative for background
245 // and positive for signal
246 // so it will be filled in the underflow bin
247 // this is very convenient for the unfolding:
248 // the unfolded result will contain the number of background
249 // events in the underflow bin
250
251 // generated MC distribution (for comparison only)
252 histMgenMC->Fill(mGen,luminosityData/luminosityMC);
253 // reconstructed MC distribution (for comparison only)
254 histMdetMC->Fill(mDet,luminosityData/luminosityMC);
255
256 // matrix describing how the generator input migrates to the
257 // reconstructed level. Unfolding input.
258 // NOTE on underflow/overflow bins:
259 // (1) the detector level under/overflow bins are used for
260 // normalisation ("efficiency" correction)
261 // in our toy example, these bins are populated from tails
262 // of the initial MC distribution.
263 // (2) the generator level underflow/overflow bins are
264 // unfolded. In this example:
265 // underflow bin: background events reconstructed in the detector
266 // overflow bin: signal events generated at masses > xmaxDet
267 // for the unfolded result these bins will be filled
268 // -> the background normalisation will be contained in the underflow bin
269 histMdetGenMC->Fill(mDet,mGen,luminosityData/luminosityMC);
270 }
271
272 //============================================
273 // generate alternative MC
274 // this will be used to derive a systematic error due to MC
275 // parameter uncertainties
276 TH2D *histMdetGenSysMC=new TH2D("MdetgenSysMC",";mass(det);mass(gen)",
277 nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
278 neventMC=rnd->Poisson(luminosityMC*crossSection);
279 for(Int_t i=0;i<neventMC;i++) {
280 Double_t mGen=GenerateEvent
281 (0.5, // relative fraction of background
282 3.6, // peak position in MC with systematic shift
283 0.15); // peak width in MC
284 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
285 histMdetGenSysMC->Fill(mDet,mGen,luminosityData/luminosityMC);
286 }
287
288 //============================================
289 // generate data distribution
290 //
291 TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen);
292 TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet);
293 Int_t neventData=rnd->Poisson(luminosityData*crossSection);
294 for(Int_t i=0;i<neventData;i++) {
295 Double_t mGen=GenerateEvent(0.4, // relative fraction of background
296 3.8, // peak position in data
297 0.15); // peak width in data
298 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
299 // generated data mass for comparison plots
300 // for real data, we do not have this histogram
301 histMgenData->Fill(mGen);
302
303 // reconstructed mass, unfolding input
304 histMdetData->Fill(mDet);
305 }
306
307 //=========================================================================
308 // divide by bin width to get density distributions
309 TH1D *histDensityGenData=new TH1D("DensityGenData",";mass(gen)",
310 nGen,xminGen,xmaxGen);
311 TH1D *histDensityGenMC=new TH1D("DensityGenMC",";mass(gen)",
312 nGen,xminGen,xmaxGen);
313 for(Int_t i=1;i<=nGen;i++) {
314 histDensityGenData->SetBinContent(i,histMgenData->GetBinContent(i)/
315 histMgenData->GetBinWidth(i));
316 histDensityGenMC->SetBinContent(i,histMgenMC->GetBinContent(i)/
317 histMgenMC->GetBinWidth(i));
318 }
319
320 //=========================================================================
321 // set up the unfolding
322 // define migration matrix
323 TUnfoldDensity unfold(histMdetGenMC,TUnfold::kHistMapOutputVert);
324
325 // define input and bias scheme
326 // do not use the bias, because MC peak may be at the wrong place
327 // watch out for error codes returned by the SetInput method
328 // errors larger or equal 10000 are fatal:
329 // the data points specified as input are not sufficient to constrain the
330 // unfolding process
331 if(unfold.SetInput(histMdetData)>=10000) {
332 std::cout<<"Unfolding result may be wrong\n";
333 }
334
335 //========================================================================
336 // the unfolding is done here
337 //
338 // scan L curve and find best point
339 Int_t nScan=30;
340 // use automatic L-curve scan: start with taumin=taumax=0.0
341 Double_t tauMin=0.0;
342 Double_t tauMax=0.0;
343 Int_t iBest;
344 TSpline *logTauX,*logTauY;
345 TGraph *lCurve;
346
347 // if required, report Info messages (for debugging the L-curve scan)
348#ifdef VERBOSE_LCURVE_SCAN
349 Int_t oldinfo=gErrorIgnoreLevel;
351#endif
352 // this method scans the parameter tau and finds the kink in the L curve
353 // finally, the unfolding is done for the best choice of tau
354 iBest=unfold.ScanLcurve(nScan,tauMin,tauMax,&lCurve,&logTauX,&logTauY);
355
356 // if required, switch to previous log-level
357#ifdef VERBOSE_LCURVE_SCAN
358 gErrorIgnoreLevel=oldinfo;
359#endif
360
361 //==========================================================================
362 // define a correlated systematic error
363 // for example, assume there is a 10% correlated error for all reconstructed
364 // masses larger than 7
365 Double_t SYS_ERROR1_MSTART=6;
366 Double_t SYS_ERROR1_SIZE=0.1;
367 TH2D *histMdetGenSys1=new TH2D("Mdetgensys1",";mass(det);mass(gen)",
368 nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
369 for(Int_t i=0;i<=nDet+1;i++) {
370 if(histMdetData->GetBinCenter(i)>=SYS_ERROR1_MSTART) {
371 for(Int_t j=0;j<=nGen+1;j++) {
372 histMdetGenSys1->SetBinContent(i,j,SYS_ERROR1_SIZE);
373 }
374 }
375 }
376 unfold.AddSysError(histMdetGenSysMC,"SYSERROR_MC",TUnfold::kHistMapOutputVert,
378 unfold.AddSysError(histMdetGenSys1,"SYSERROR1",TUnfold::kHistMapOutputVert,
380
381 //==========================================================================
382 // print some results
383 //
384 std::cout<<"tau="<<unfold.GetTau()<<"\n";
385 std::cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
386 <<" / "<<unfold.GetNdf()<<"\n";
387 std::cout<<"chi**2(sys)="<<unfold.GetChi2Sys()<<"\n";
388
389
390 //==========================================================================
391 // create graphs with one point to visualize the best choice of tau
392 //
393 Double_t t[1],x[1],y[1];
394 logTauX->GetKnot(iBest,t[0],x[0]);
395 logTauY->GetKnot(iBest,t[0],y[0]);
396 TGraph *bestLcurve=new TGraph(1,x,y);
397 TGraph *bestLogTauLogChi2=new TGraph(1,t,x);
398
399 //==========================================================================
400 // retrieve results into histograms
401
402 // get unfolded distribution
403 TH1 *histMunfold=unfold.GetOutput("Unfolded");
404
405 // get unfolding result, folded back
406 TH1 *histMdetFold=unfold.GetFoldedOutput("FoldedBack");
407
408 // get error matrix (input distribution [stat] errors only)
409 // TH2D *histEmatData=unfold.GetEmatrix("EmatData");
410
411 // get total error matrix:
412 // migration matrix uncorrelated and correlated systematic errors
413 // added in quadrature to the data statistical errors
414 TH2 *histEmatTotal=unfold.GetEmatrixTotal("EmatTotal");
415
416 // create data histogram with the total errors
417 TH1D *histTotalError=
418 new TH1D("TotalError",";mass(gen)",nGen,xminGen,xmaxGen);
419 for(Int_t bin=1;bin<=nGen;bin++) {
420 histTotalError->SetBinContent(bin,histMunfold->GetBinContent(bin));
421 histTotalError->SetBinError
422 (bin,TMath::Sqrt(histEmatTotal->GetBinContent(bin,bin)));
423 }
424
425 // get global correlation coefficients
426 // for this calculation one has to specify whether the
427 // underflow/overflow bins are included or not
428 // default: include all bins
429 // here: exclude underflow and overflow bins
430 TH1 *histRhoi=unfold.GetRhoItotal("rho_I",
431 nullptr, // use default title
432 nullptr, // all distributions
433 "*[UO]", // discard underflow and overflow bins on all axes
434 kTRUE, // use original binning
435 &gHistInvEMatrix // store inverse of error matrix
436 );
437
438 //======================================================================
439 // fit Breit-Wigner shape to unfolded data, using the full error matrix
440 // here we use a "user" chi**2 function to take into account
441 // the full covariance matrix
442
443 gFitter=TVirtualFitter::Fitter(histMunfold);
444 gFitter->SetFCN(chisquare_corr);
445
446 TF1 *bw=new TF1("bw",bw_func,xminGen,xmaxGen,3);
447 bw->SetParameter(0,1000.);
448 bw->SetParameter(1,3.8);
449 bw->SetParameter(2,0.2);
450
451 // for (wrong!) fitting without correlations, drop the option "U"
452 // here.
453 histMunfold->Fit(bw,"UE");
454
455 //=====================================================================
456 // plot some histograms
458 output.Divide(3,2);
459
460 // Show the matrix which connects input and output
461 // There are overflow bins at the bottom, not shown in the plot
462 // These contain the background shape.
463 // The overflow bins to the left and right contain
464 // events which are not reconstructed. These are necessary for proper MC
465 // normalisation
466 output.cd(1);
467 histMdetGenMC->Draw("BOX");
468
469 // draw generator-level distribution:
470 // data (red) [for real data this is not available]
471 // MC input (black) [with completely wrong peak position and shape]
472 // unfolded data (blue)
473 output.cd(2);
474 histTotalError->SetLineColor(kBlue);
475 histTotalError->Draw("E");
476 histMunfold->SetLineColor(kGreen);
477 histMunfold->Draw("SAME E1");
478 histDensityGenData->SetLineColor(kRed);
479 histDensityGenData->Draw("SAME");
480 histDensityGenMC->Draw("SAME HIST");
481
482 // show detector level distributions
483 // data (red)
484 // MC (black) [with completely wrong peak position and shape]
485 // unfolded data (blue)
486 output.cd(3);
487 histMdetFold->SetLineColor(kBlue);
488 histMdetFold->Draw();
489 histMdetMC->Draw("SAME HIST");
490
491 TH1 *histInput=unfold.GetInput("Minput",";mass(det)");
492
493 histInput->SetLineColor(kRed);
494 histInput->Draw("SAME");
495
496 // show correlation coefficients
497 output.cd(4);
498 histRhoi->Draw();
499
500 // show tau as a function of chi**2
501 output.cd(5);
502 logTauX->Draw();
503 bestLogTauLogChi2->SetMarkerColor(kRed);
504 bestLogTauLogChi2->Draw("*");
505
506 // show the L curve
507 output.cd(6);
508 lCurve->Draw("AL");
509 bestLcurve->SetMarkerColor(kRed);
510 bestLcurve->Draw("*");
511
512 output.SaveAs("testUnfold1.ps");
513
514 return 0;
515}
#define f(i)
Definition RSha256.hxx:104
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
@ kRed
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
constexpr Int_t kInfo
Definition TError.h:45
Int_t gErrorIgnoreLevel
Error handling routines.
Definition TError.cxx:31
R__EXTERN TStyle * gStyle
Definition TStyle.h:436
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:23
1-Dim function class
Definition TF1.h:233
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:3683
virtual Int_t GetNpar() const
Definition TF1.h:509
virtual void SetNumberFitPoints(Int_t npfits)
Definition TF1.h:652
virtual void InitArgs(const Double_t *x, const Double_t *params)
Initialize parameters addresses.
Definition TF1.cxx:2482
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr)
Evaluate function with given coordinates and parameters.
Definition TF1.cxx:1468
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition TF1.cxx:3692
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:667
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:833
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:671
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
Definition TH1.cxx:9174
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:3906
virtual Int_t GetNbinsX() const
Definition TH1.h:298
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:9239
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3346
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3068
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:6732
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:9255
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5090
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width for 1D histogram.
Definition TH1.cxx:9196
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:358
Service class for 2-D histogram classes.
Definition TH2.h:30
void SetBinContent(Int_t bin, Double_t content) override
Set bin content.
Definition TH2.cxx:2613
Double_t GetBinContent(Int_t binx, Int_t biny) const override
Definition TH2.h:94
Int_t Fill(Double_t) override
Invalid Fill method.
Definition TH2.cxx:393
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Automatic pad generation by division.
Definition TPad.cxx:1249
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:275
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom.cxx:559
virtual ULong64_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
Definition TRandom.cxx:404
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
void Draw(Option_t *option="") override
Draw this function with its current attributes.
Definition TSpline.cxx:101
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:1593
An algorithm to unfold distributions from detector to truth level.
@ kSysErrModeRelative
matrix gives the relative shifts
Definition TUnfoldSys.h:112
@ kSysErrModeMatrix
matrix is an alternative to the default matrix, the errors are the difference to the original matrix
Definition TUnfoldSys.h:108
@ kHistMapOutputVert
truth level on y-axis of the response matrix
Definition TUnfold.h:149
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)
Returns the square root of x.
Definition TMath.h:666
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:725
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Tan(Double_t)
Returns the tangent of an angle of x radians.
Definition TMath.h:604
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
static void output()