Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
testUnfold3.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_unfold
3/// \notebook
4/// Simple Test program for the class TUnfoldDensity.
5///
6/// 1-dimensional unfolding with background subtraction
7///
8/// the goal is to unfold the underlying "true" distribution of a variable Pt
9///
10/// the reconstructed Pt is measured in 24 bins from 4 to 28
11/// the generator-level Pt is unfolded into 10 bins from 6 to 26
12/// - plus underflow bin from 0 to 6
13/// - plus overflow bin above 26
14/// there are two background sources bgr1 and bgr2
15/// the signal has a finite trigger efficiency at a threshold of 8 GeV
16///
17/// one type of systematic error is studied, where the signal parameters are
18/// changed
19///
20/// Finally, the unfolding is compared to a "bin-by-bin" correction method
21///
22/// \macro_output
23/// \macro_code
24///
25/// **Version 17.6, in parallel to changes in TUnfold**
26///
27/// #### History:
28/// - Version 17.5, in parallel to changes in TUnfold
29/// - Version 17.4, in parallel to changes in TUnfold
30/// - Version 17.3, in parallel to changes in TUnfold
31/// - Version 17.2, in parallel to changes in TUnfold
32/// - Version 17.1, in parallel to changes in TUnfold
33/// - Version 17.0, change to use the class TUnfoldDensity
34/// - Version 16.1, parallel to changes in TUnfold
35/// - Version 16.0, parallel to changes in TUnfold
36/// - Version 15, simple example including background subtraction
37///
38/// This file is part of TUnfold.
39///
40/// TUnfold is free software: you can redistribute it and/or modify
41/// it under the terms of the GNU General Public License as published by
42/// the Free Software Foundation, either version 3 of the License, or
43/// (at your option) any later version.
44///
45/// TUnfold is distributed in the hope that it will be useful,
46/// but WITHOUT ANY WARRANTY; without even the implied warranty of
47/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
48/// GNU General Public License for more details.
49///
50/// You should have received a copy of the GNU General Public License
51/// along with TUnfold. If not, see <http://www.gnu.org/licenses/>.
52///
53/// \author Stefan Schmitt DESY, 14.10.2008
54
55#include <TMath.h>
56#include <TCanvas.h>
57#include <TRandom3.h>
58#include <TF1.h>
59#include <TStyle.h>
60#include <TVector.h>
61#include <TGraph.h>
62
63#include "TUnfoldDensity.h"
64
65TRandom *rnd=nullptr;
66
67Double_t GenerateEvent(const Double_t *parm,
68 const Double_t *triggerParm,
72{
73 // generate an event
74 // input:
75 // parameters for the event generator
76 // return value:
77 // reconstructed Pt
78 // output to pointers:
79 // integrated luminosity
80 // several variables only accessible on generator level
81 //
82 // the parm array defines the physical parameters
83 // parm[0]: background source 1 fraction
84 // parm[1]: background source 2 fraction
85 // parm[2]: lower limit of generated Pt distribution
86 // parm[3]: upper limit of generated Pt distribution
87 // parm[4]: exponent for Pt distribution signal
88 // parm[5]: exponent for Pt distribution background 1
89 // parm[6]: exponent for Pt distribution background 2
90 // parm[7]: resolution parameter a goes with sqrt(Pt)
91 // parm[8]: resolution parameter b goes with Pt
92 // triggerParm[0]: trigger threshold turn-on position
93 // triggerParm[1]: trigger threshold turn-on width
94 // triggerParm[2]: trigger efficiency for high pt
95 //
96 // intLumi is advanced bu 1 for each *generated* event
97 // for data, several events may be generated, until one passes the trigger
98 //
99 // some generator-level quantities are also returned:
100 // triggerFlag: whether the event passed the trigger threshold
101 // ptGen: the generated pt
102 // iType: which type of process was simulated
103 //
104 // the "triggerFlag" also has another meaning:
105 // if(triggerFlag==0) only triggered events are returned
106 //
107 // Usage to generate data events
108 // ptObs=GenerateEvent(parm,triggerParm,0,0,0)
109 //
110 // Usage to generate MC events
111 // ptGen=GenerateEvent(parm,triggerParm,&triggerFlag,&ptGen,&iType);
112 //
115 do {
116 Int_t itype;
118 Double_t f=rnd->Rndm();
119 // decide whether this is background or signal
120 itype=0;
121 if(f<parm[0]) itype=1;
122 else if(f<parm[0]+parm[1]) itype=2;
123 // generate Pt according to distribution pow(ptgen,a)
124 // get exponent
125 Double_t a=parm[4+itype];
126 Double_t a1=a+1.0;
127 Double_t t=rnd->Rndm();
128 if(a1 == 0.0) {
129 Double_t x0=TMath::Log(parm[2]);
130 ptgen=TMath::Exp(t*(TMath::Log(parm[3])-x0)+x0);
131 } else {
132 Double_t x0=pow(parm[2],a1);
133 ptgen=pow(t*(pow(parm[3],a1)-x0)+x0,1./a1);
134 }
135 if(iType) *iType=itype;
136 if(ptGen) *ptGen=ptgen;
137
138 // smearing in Pt with large asymmetric tail
141 ptObs=rnd->BreitWigner(ptgen,sigma);
142
143 // decide whether event was triggered
146 isTriggered= rnd->Rndm()<triggerProb;
147 (*intLumi) ++;
148 } while((!triggerFlag) && (!isTriggered));
149 // return trigger decision
151 return ptObs;
152}
153
154//int main(int argc, char *argv[])
155void testUnfold3()
156{
157 // switch on histogram errors
159
160 // show fit result
161 gStyle->SetOptFit(1111);
162
163 // random generator
164 rnd=new TRandom3();
165
166 // data and MC luminosities
167 Double_t const lumiData= 30000;
168 Double_t const lumiMC =1000000;
169
170 //========================
171 // Step 1: define binning, book histograms
172
173 // reconstructed pt (fine binning)
174 Int_t const nDet=24;
175 Double_t const xminDet=4.0;
176 Double_t const xmaxDet=28.0;
177
178 // generated pt (coarse binning)
179 Int_t const nGen=10;
180 Double_t const xminGen= 6.0;
181 Double_t const xmaxGen=26.0;
182
183 //==================================================================
184 // book histograms
185 // (1) unfolding input: binning scheme is fine for detector,coarse for gen
186 // histUnfoldInput : reconstructed data, binning for unfolding
187 // histUnfoldMatrix : generated vs reconstructed distribution
188 // histUnfoldBgr1 : background source1, as predicted by MC
189 // histUnfoldBgr2 : background source2, as predicted by MC
190 // for systematic studies
191 // histUnfoldMatrixSys : generated vs reconstructed with different signal
192 //
193 // (2) histograms required for bin-by-bin method
194 // histDetDATAbbb : reconstructed data for bin-by-bin
195 // histDetMCbbb : reconstructed MC for bin-by-bin
196 // histDetMCBgr1bbb : reconstructed MC bgr1 for bin-by-bin
197 // histDetMCBgr2bbb : reconstructed MC bgr2 for bin-by-bin
198 // histDetMCBgrPtbbb : reconstructed MC bgr from low/high pt for bin-by-bin
199 // histGenMCbbb : generated MC truth
200 // for systematic studies
201 // histDetMCSysbbb : reconstructed with changed trigger
202 // histGenMCSysbbb : generated MC truth
203 //
204 // (3) monitoring and control
205 // histGenData : data truth for bias tests
206 // histDetMC : MC prediction
207
208 // (1) create histograms required for unfolding
210 new TH1D("unfolding input rec",";ptrec",nDet,xminDet,xmaxDet);
212 new TH2D("unfolding matrix",";ptgen;ptrec",
215 new TH1D("unfolding bgr1 rec",";ptrec",nDet,xminDet,xmaxDet);
217 new TH1D("unfolding bgr2 rec",";ptrec",nDet,xminDet,xmaxDet);
219 new TH2D("unfolding matrix sys",";ptgen;ptrec",
221
222 // (2) histograms required for bin-by-bin
224 new TH1D("bbb input rec",";ptrec",nGen,xminGen,xmaxGen);
226 new TH1D("bbb signal rec",";ptrec",nGen,xminGen,xmaxGen);
228 new TH1D("bbb bgr1 rec",";ptrec",nGen,xminGen,xmaxGen);
230 new TH1D("bbb bgr2 rec",";ptrec",nGen,xminGen,xmaxGen);
232 new TH1D("bbb bgrPt rec",";ptrec",nGen,xminGen,xmaxGen);
234 new TH1D("bbb signal gen",";ptgen",nGen,xminGen,xmaxGen);
236 new TH1D("bbb signal sys rec",";ptrec",nGen,xminGen,xmaxGen);
238 new TH1D("bbb bgrPt sys rec",";ptrec",nGen,xminGen,xmaxGen);
240 new TH1D("bbb signal sys gen",";ptgen",nGen,xminGen,xmaxGen);
241
242 // (3) control histograms
244 new TH1D("DATA truth gen",";ptgen",nGen,xminGen,xmaxGen);
246 new TH1D("MC prediction rec",";ptrec",nDet,xminDet,xmaxDet);
247 // ==============================================================
248 // Step 2: generate data distributions
249 //
250 // data parameters: in real life these are unknown
251 static Double_t parm_DATA[]={
252 0.05, // fraction of background 1 (on generator level)
253 0.05, // fraction of background 2 (on generator level)
254 3.5, // lower Pt cut (generator level)
255 100.,// upper Pt cut (generator level)
256 -3.0,// signal exponent
257 0.1, // background 1 exponent
258 -0.5, // background 2 exponent
259 0.2, // energy resolution a term
260 0.01, // energy resolution b term
261 };
262 static Double_t triggerParm_DATA[]={8.0, // threshold is 8 GeV
263 4.0, // width is 4 GeV
264 0.95 // high Pt efficiency os 95%
265 };
266
267 Double_t intLumi=0.0;
268 while(intLumi<lumiData) {
274 if(isTriggered) {
275 // (1) histogram for unfolding
276 histUnfoldInput->Fill(ptObs);
277
278 // (2) histogram for bin-by-bin
279 histBbbInput->Fill(ptObs);
280 }
281 // (3) monitoring
282 if(iTypeGen==0) histDataTruth->Fill(ptGen);
283 }
284
285 // ==============================================================
286 // Step 3: generate default MC distributions
287 //
288 // MC parameters
289 // default settings
290 static Double_t parm_MC[]={
291 0.05, // fraction of background 1 (on generator level)
292 0.05, // fraction of background 2 (on generator level)
293 3.5, // lower Pt cut (generator level)
294 100.,// upper Pt cut (generator level)
295 -4.0,// signal exponent !!! steeper than in data
296 // to illustrate bin-by-bin bias
297 0.1, // background 1 exponent
298 -0.5, // background 2 exponent
299 0.2, // energy resolution a term
300 0.01, // energy resolution b term
301 };
302 static Double_t triggerParm_MC[]={8.0, // threshold is 8 GeV
303 4.0, // width is 4 GeV
304 0.95 // high Pt efficiency is 95%
305 };
306
307 // weighting factor to accomodate for the different luminosity in data and MC
309 intLumi=0.0;
310 while(intLumi<lumiMC) {
315 &ptGen,&iTypeGen);
316 if(!isTriggered) ptObs=0.0;
317
318 // (1) distribution required for unfolding
319
320 if(iTypeGen==0) {
322 } else if(iTypeGen==1) {
324 } else if(iTypeGen==2) {
326 }
327
328 // (2) distributions for bbb
329 if(iTypeGen==0) {
330 if((ptGen>=xminGen)&&(ptGen<xmaxGen)) {
332 } else {
334 }
336 } else if(iTypeGen==1) {
338 } else if(iTypeGen==2) {
340 }
341
342 // (3) control distribution
343 histDetMC ->Fill(ptObs,lumiWeight);
344 }
345
346 // ==============================================================
347 // Step 4: generate MC distributions for systematic study
348 // example: study dependence on initial signal shape
349 // -> BGR distributions do not change
350 static Double_t parm_MC_SYS[]={
351 0.05, // fraction of background: unchanged
352 0.05, // fraction of background: unchanged
353 3.5, // lower Pt cut (generator level)
354 100.,// upper Pt cut (generator level)
355 -2.0, // signal exponent changed
356 0.1, // background 1 exponent
357 -0.5, // background 2 exponent
358 0.2, // energy resolution a term
359 0.01, // energy resolution b term
360 };
361
362 intLumi=0.0;
363 while(intLumi<lumiMC) {
369 if(!isTriggered) ptObs=0.0;
370
371 // (1) for unfolding
372 if(iTypeGen==0) {
374 }
375
376 // (2) for bin-by-bin
377 if(iTypeGen==0) {
378 if((ptGen>=xminGen)&&(ptGen<xmaxGen)) {
380 } else {
381 histBbbBgrPtSys->Fill(ptObs);
382 }
384 }
385 }
386
387 // this method is new in version 16 of TUnfold
388 cout<<"TUnfold version is "<<TUnfold::GetTUnfoldVersion()<<"\n";
389
390 //========================
391 // Step 5: unfolding
392 //
393
394 // here one has to decide about the regularisation scheme
395
396 // regularize curvature
399
400 // preserve the area
403
404 // bin content is divided by the bin width
407
408 // set up matrix of migrations
411
412 // define the input vector (the measured data distribution)
413 unfold.SetInput(histUnfoldInput);
414
415 // subtract background, normalized to data luminosity
416 // with 10% scale error each
419 unfold.SubtractBackground(histUnfoldBgr1,"background1",scale_bgr,dscale_bgr);
420 unfold.SubtractBackground(histUnfoldBgr2,"background2",scale_bgr,dscale_bgr);
421
422 // add systematic error
423 unfold.AddSysError(histUnfoldMatrixSys,"signalshape_SYS",
426
427 // run the unfolding
428 Int_t nScan=30;
430 TGraph *lCurve;
431 // this method scans the parameter tau and finds the kink in the L curve
432 // finally, the unfolding is done for the best choice of tau
433 Int_t iBest=unfold.ScanLcurve(nScan,0.,0.,&lCurve,&logTauX,&logTauY);
434
435 cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
436 <<" / "<<unfold.GetNdf()<<"\n";
437
438 // save graphs with one point to visualize best choice of tau
439 Double_t t[1],x[1],y[1];
440 logTauX->GetKnot(iBest,t[0],x[0]);
441 logTauY->GetKnot(iBest,t[0],y[0]);
442 TGraph *bestLcurve=new TGraph(1,x,y);
444
445
446 //===========================
447 // Step 6: retrieve unfolding results
448
449 // get unfolding output
450 // includes the statistical and background errors
451 // but not the other systematic uncertainties
452 TH1 *histUnfoldOutput=unfold.GetOutput("PT(unfold,stat+bgrerr)");
453
454 // retrieve error matrix of statistical errors
455 TH2 *histEmatStat=unfold.GetEmatrixInput("unfolding stat error matrix");
456
457 // retrieve full error matrix
458 // This includes all systematic errors
459 TH2 *histEmatTotal=unfold.GetEmatrixTotal("unfolding total error matrix");
460
461 // create two copies of the unfolded data, one with statistical errors
462 // the other with total errors
463 TH1 *histUnfoldStat=new TH1D("PT(unfold,staterr)",";Pt(gen)",
465 TH1 *histUnfoldTotal=new TH1D("PT(unfold,totalerr)",";Pt(gen)",
467 for(Int_t i=0;i<nGen+2;i++) {
468 Double_t c=histUnfoldOutput->GetBinContent(i);
469
470 // histogram with unfolded data and stat errors
471 histUnfoldStat->SetBinContent(i,c);
472 histUnfoldStat->SetBinError
473 (i,TMath::Sqrt(histEmatStat->GetBinContent(i,i)));
474
475 // histogram with unfolded data and total errors
476 histUnfoldTotal->SetBinContent(i,c);
477 histUnfoldTotal->SetBinError
478 (i,TMath::Sqrt(histEmatTotal->GetBinContent(i,i)));
479 }
480
481 // create histogram with correlation matrix
482 TH2D *histCorr=new TH2D("Corr(total)",";Pt(gen);Pt(gen)",
484 for(Int_t i=0;i<nGen+2;i++) {
485 Double_t ei,ej;
486 ei=TMath::Sqrt(histEmatTotal->GetBinContent(i,i));
487 if(ei<=0.0) continue;
488 for(Int_t j=0;j<nGen+2;j++) {
489 ej=TMath::Sqrt(histEmatTotal->GetBinContent(j,j));
490 if(ej<=0.0) continue;
491 histCorr->SetBinContent(i,j,histEmatTotal->GetBinContent(i,j)/ei/ej);
492 }
493 }
494
495 // retrieve bgr source 1
496 TH1 *histDetNormBgr1=unfold.GetBackground("bgr1 normalized",
497 "background1");
498 TH1 *histDetNormBgrTotal=unfold.GetBackground("bgr total normalized");
499
500 //========================
501 // Step 7: plots
502
504 output.Divide(3,2);
505 output.cd(1);
506 // data, MC prediction, background
507 histUnfoldInput->SetMinimum(0.0);
508 histUnfoldInput->Draw("E");
509 histDetMC->SetMinimum(0.0);
510 histDetMC->SetLineColor(kBlue);
511 histDetNormBgrTotal->SetLineColor(kRed);
512 histDetNormBgr1->SetLineColor(kCyan);
513 histDetMC->Draw("SAME HIST");
514 histDetNormBgr1->Draw("SAME HIST");
515 histDetNormBgrTotal->Draw("SAME HIST");
516
517 output.cd(2);
518 // unfolded data, data truth, MC truth
519 histUnfoldTotal->SetMinimum(0.0);
520 histUnfoldTotal->SetMaximum(histUnfoldTotal->GetMaximum()*1.5);
521 // outer error: total error
522 histUnfoldTotal->Draw("E");
523 // middle error: stat+bgr
524 histUnfoldOutput->Draw("SAME E1");
525 // inner error: stat only
526 histUnfoldStat->Draw("SAME E1");
527 histDataTruth->Draw("SAME HIST");
528 histBbbSignalGen->SetLineColor(kBlue);
529 histBbbSignalGen->Draw("SAME HIST");
530
531 output.cd(3);
532 // unfolding matrix
533 histUnfoldMatrix->SetLineColor(kBlue);
534 histUnfoldMatrix->Draw("BOX");
535
536 // show tau as a function of chi**2
537 output.cd(4);
538 logTauX->Draw();
539 bestLogTauLogChi2->SetMarkerColor(kRed);
540 bestLogTauLogChi2->Draw("*");
541
542 // show the L curve
543 output.cd(5);
544 lCurve->Draw("AL");
545 bestLcurve->SetMarkerColor(kRed);
546 bestLcurve->Draw("*");
547
548 // show correlation matrix
549 output.cd(6);
550 histCorr->Draw("BOX");
551
552 output.SaveAs("testUnfold3.ps");
553
554 //============================================================
555 // step 8: compare results to the so-called bin-by-bin "correction"
556
557 std::cout<<"bin truth result (stat) (bgr) (sys)\n";
558 std::cout<<"===+=====+=========+========+========+=======\n";
559 for(Int_t i=1;i<=nGen;i++) {
560 // data contribution in this bin
561 Double_t data=histBbbInput->GetBinContent(i);
562 Double_t errData=histBbbInput->GetBinError(i);
563
564 // subtract background contributions
566 - scale_bgr*(histBbbBgr1->GetBinContent(i)
567 + histBbbBgr2->GetBinContent(i)
568 + histBbbBgrPt->GetBinContent(i));
571 TMath::Power(dscale_bgr*histBbbBgr1->GetBinContent(i),2)+
572 TMath::Power(scale_bgr*histBbbBgr1->GetBinError(i),2)+
573 TMath::Power(dscale_bgr*histBbbBgr2->GetBinContent(i),2)+
574 TMath::Power(scale_bgr*histBbbBgr2->GetBinError(i),2)+
575 TMath::Power(dscale_bgr*histBbbBgrPt->GetBinContent(i),2)+
576 TMath::Power(scale_bgr*histBbbBgrPt->GetBinError(i),2));
577 // "correct" the data, using the Monte Carlo and neglecting off-diagonals
578 Double_t fCorr=(histBbbSignalGen->GetBinContent(i)/
579 histBbbSignalRec->GetBinContent(i));
581 // stat only error
583 // stat plus background subtraction
585
586 // estimate systematic error by repeating the exercise
587 // using the MC with systematic shifts applied
588 Double_t fCorr_sys=(histBbbSignalGenSys->GetBinContent(i)/
589 histBbbSignalRecSys->GetBinContent(i));
591
592 // add systematic shift quadratically and get total error
596
597 // get results from real unfolding
598 Double_t data_unfold= histUnfoldStat->GetBinContent(i);
602
603 // compare
604
605 std::cout<<TString::Format
606 ("%3d %5.0f %8.1f +/-%5.1f +/-%5.1f +/-%5.1f (unfolding)",
607 i,histDataTruth->GetBinContent(i),data_unfold,
616 std::cout<<TString::Format
617 (" %8.1f +/-%5.1f +/-%5.1f +/-%5.1f (bin-by-bin)",
626 <<"\n\n";
627 }
628}
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
bool Bool_t
Boolean (0=false, 1=true) (bool)
Definition RtypesCore.h:77
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
@ kRed
Definition Rtypes.h:67
@ kCyan
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
The Canvas class.
Definition TCanvas.h:23
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:926
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
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:6751
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:400
Service class for 2-D histogram classes.
Definition TH2.h:30
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
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2384
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:1594
An algorithm to unfold distributions from detector to truth level.
EDensityMode
choice of regularisation scale factors to cinstruct the matrix L
@ kDensityModeBinWidth
scale factors from multidimensional bin width
@ kSysErrModeMatrix
matrix is an alternative to the default matrix, the errors are the difference to the original matrix
Definition TUnfoldSys.h:108
static const char * GetTUnfoldVersion(void)
return a string describing the TUnfold version
Definition TUnfold.cxx:3717
EConstraint
type of extra constraint
Definition TUnfold.h:113
@ kEConstraintArea
enforce preservation of the area
Definition TUnfold.h:119
ERegMode
choice of regularisation scheme
Definition TUnfold.h:123
@ kRegModeCurvature
regularize the 2nd derivative of the output distribution
Definition TUnfold.h:135
@ kHistMapOutputHoriz
truth level on x-axis of the response matrix
Definition TUnfold.h:146
const Double_t sigma
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
Definition TMath.h:720
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:767
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:732
static void output()