Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
testUnfold7c.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Test program for the classes TUnfoldDensity and TUnfoldBinning.

A toy test of the TUnfold package

This example is documented in conference proceedings:

arXiv:1611.01927 12th Conference on Quark Confinement and the Hadron Spectrum (Confinement XII)

This is an example of unfolding a one-dimensional distribution. It compares various unfolding methods:

      matrix inversion, template fit, L-curve scan,
      iterative unfolding, etc

Further details can be found in talk by S.Schmitt at:

XII Quark Confinement and the Hadron Spectrum 29.8. - 3.9.2016 Thessaloniki, Greece statictics session (+proceedings)

The example comprises several macros

  • testUnfold7a.C create root files with TTree objects for signal, background and data
    • write files testUnfold7_signal.root testUnfold7_background.root testUnfold7_data.root
  • testUnfold7b.C loop over trees and fill histograms based on the TUnfoldBinning objects
    • read testUnfold7binning.xml testUnfold7_signal.root testUnfold7_background.root testUnfold7_data.root
    • write testUnfold7_histograms.root
  • testUnfold7c.C run the unfolding
    • read testUnfold7_histograms.root
    • write testUnfold7_result.root
    • write many histograms, to compare various unfolding methods
maximum number of iterations: 1000
0
100
200
300
400
500
600
700
800
900
1000
histOutputCtau0_binwU
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 8.93341
NDf = 0
Edm = 2.94432e-07
NCalls = 119
p0 = 2313.03 +/- 81.4707 -80.0337 +82.9445 (Minos)
p1 = 6.07588 +/- 0.0843195 -0.0839209 +0.0847407 (Minos)
p2 = 1.76119 +/- 0.0622688 -0.0617262 +0.0628078 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2313.03 6.07588 1.76119]
histOutputFtau0_binwU
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 8.33089
NDf = 0
Edm = 6.50812e-08
NCalls = 119
p0 = 2292.59 +/- 77.1869 -75.9067 +78.5004 (Minos)
p1 = 6.0667 +/- 0.0821746 -0.0817803 +0.0825885 (Minos)
p2 = 1.77123 +/- 0.0599625 -0.0594568 +0.0604685 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2292.59 6.0667 1.77123]
histOutputFAtau0_binwU
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 8.36913
NDf = 0
Edm = 4.1629e-08
NCalls = 119
p0 = 2300.96 +/- 77.0613 -75.7926 +78.3631 (Minos)
p1 = 6.06894 +/- 0.0818049 -0.0814136 +0.0822153 (Minos)
p2 = 1.77403 +/- 0.0597178 -0.0592155 +0.0602208 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2300.96 6.06894 1.77403]
histOutputFALCurve_binwU
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 22.1459
NDf = 0
Edm = 4.72718e-08
NCalls = 123
p0 = 2175.69 +/- 63.3586 -62.8717 +63.8248 (Minos)
p1 = 6.01392 +/- 0.0830555 -0.08261 +0.08352 (Minos)
p2 = 1.88219 +/- 0.0523609 -0.0516227 +0.0531036 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2175.69 6.01392 1.88219]
histOutputFArho_binwU
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 36.4398
NDf = 0
Edm = 8.28607e-06
NCalls = 122
p0 = 2088.52 +/- 54.3348 -54.0383 +54.6064 (Minos)
p1 = 5.94814 +/- 0.0813071 -0.0810537 +0.0815775 (Minos)
p2 = 1.96606 +/- 0.0470327 -0.0463947 +0.0476727 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2088.52 5.94814 1.96606]
histOutputBBBO
bad global correlation 8 -2.22045e-16
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 19.0585
NDf = 0
Edm = 1.0184e-08
NCalls = 102
p0 = 1963.05 +/- 51.5682 -51.2648 +51.8646 (Minos)
p1 = 5.90214 +/- 0.0816763 -0.0809708 +0.0824027 (Minos)
p2 = 2.07111 +/- 0.0460756 -0.045646 +0.0465268 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1963.05 5.90214 2.07111]
histOutputAgorep0
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 197.953
NDf = 0
Edm = 5.71447e-06
NCalls = 101
p0 = 1707.06 +/- 30.628 -30.6773 +30.5555 (Minos)
p1 = 5.62399 +/- 0.0532142 -0.052652 +0.0537917 (Minos)
p2 = 2.24364 +/- 0.020155 -0.0199327 +0.0203847 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1707.06 5.62399 2.24364]
histOutputAgorep1
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 118.197
NDf = 0
Edm = 7.06309e-08
NCalls = 120
p0 = 1766.67 +/- 33.976 -33.9516 +33.9857 (Minos)
p1 = 5.8135 +/- 0.0699476 -0.0693366 +0.0706199 (Minos)
p2 = 2.27527 +/- 0.0262139 -0.0259602 +0.0264935 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1766.67 5.8135 2.27527]
histOutputAgorep2
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 103.621
NDf = 0
Edm = 1.74226e-07
NCalls = 123
p0 = 1794.2 +/- 35.0632 -35.0418 +35.0726 (Minos)
p1 = 5.89902 +/- 0.0771589 -0.0765447 +0.077804 (Minos)
p2 = 2.27417 +/- 0.0289737 -0.028679 +0.0292851 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1794.2 5.89902 2.27417]
histOutputAgorep3
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 96.8081
NDf = 0
Edm = 2.93973e-06
NCalls = 119
p0 = 1816.25 +/- 35.9125 -35.9098 +35.9051 (Minos)
p1 = 5.93729 +/- 0.0805321 -0.0799526 +0.0811401 (Minos)
p2 = 2.25889 +/- 0.0305052 -0.0301629 +0.0308666 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1816.25 5.93729 2.25889]
histOutputAgorep4
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 91.6257
NDf = 0
Edm = 5.08799e-07
NCalls = 123
p0 = 1836.74 +/- 36.7459 -36.7349 +36.7473 (Minos)
p1 = 5.9568 +/- 0.0823367 -0.0816809 +0.0830262 (Minos)
p2 = 2.2385 +/- 0.0316079 -0.0312164 +0.0320162 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1836.74 5.9568 2.2385]
histOutputAgorep5
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 86.8941
NDf = 0
Edm = 5.19992e-07
NCalls = 120
p0 = 1856.65 +/- 37.6417 -37.6077 +37.665 (Minos)
p1 = 5.96958 +/- 0.0835144 -0.082935 +0.0841269 (Minos)
p2 = 2.2164 +/- 0.0326294 -0.0322428 +0.0330273 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1856.65 5.96958 2.2164]
histOutputAgorep6
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 82.3357
NDf = 0
Edm = 3.22307e-07
NCalls = 115
p0 = 1876.57 +/- 38.6709 -38.6446 +38.6871 (Minos)
p1 = 5.98002 +/- 0.0844752 -0.0837658 +0.0852298 (Minos)
p2 = 2.1935 +/- 0.0337655 -0.0333551 +0.0341892 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1876.57 5.98002 2.1935]
histOutputAgorep7
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 77.8736
NDf = 0
Edm = 2.18856e-07
NCalls = 116
p0 = 1896.9 +/- 39.939 -39.8492 +40.0236 (Minos)
p1 = 5.99052 +/- 0.0854688 -0.0847383 +0.0862587 (Minos)
p2 = 2.16994 +/- 0.0352059 -0.0348797 +0.035563 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1896.9 5.99052 2.16994]
histOutputAgorep8
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 73.4731
NDf = 0
Edm = 7.26758e-08
NCalls = 103
p0 = 1918.35 +/- 41.6429 -41.4587 +41.8498 (Minos)
p1 = 6.0021 +/- 0.0867252 -0.0858723 +0.0876869 (Minos)
p2 = 2.14524 +/- 0.0372357 -0.0371887 +0.0374098 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1918.35 6.0021 2.14524]
histOutputAgorep9
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 69.0987
NDf = 0
Edm = 5.25891e-08
NCalls = 140
p0 = 1941.89 +/- 44.2607 -43.7871 +44.9422 (Minos)
p1 = 6.01603 +/- 0.0887007 -0.0875698 +0.0901215 (Minos)
p2 = 2.11836 +/- 0.0405096 -0.0414649 +0.0401739 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 1941.89 6.01603 2.11836]
histOutputAgorep10
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 42.5138
NDf = 0
Edm = 5.85033e-07
NCalls = 131
p0 = 2382.79 +/- 63.0985 -63.4333 +62.8409 (Minos)
p1 = 6.23206 +/- 0.0749671 -0.0743234 +0.07562 (Minos)
p2 = 1.6911 +/- 0.0410841 -0.0398794 +0.0425078 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2382.79 6.23206 1.6911]
histOutputAgorep11
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 38.7947
NDf = 0
Edm = 3.1702e-08
NCalls = 115
p0 = 2389.25 +/- 65.1256 -65.6257 +64.7417 (Minos)
p1 = 6.21773 +/- 0.0753106 -0.0746172 +0.0760137 (Minos)
p2 = 1.6883 +/- 0.0428048 -0.0413996 +0.0444921 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2389.25 6.21773 1.6883]
histOutputAgorep12
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 35.7239
NDf = 0
Edm = 5.37905e-07
NCalls = 98
p0 = 2392.49 +/- 67.1401 -67.7426 +66.6974 (Minos)
p1 = 6.20537 +/- 0.0756438 -0.0750129 +0.0762805 (Minos)
p2 = 1.68766 +/- 0.0446165 -0.0430544 +0.0465276 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2392.49 6.20537 1.68766]
histOutputAgorep13
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 33.1305
NDf = 0
Edm = 4.98975e-09
NCalls = 108
p0 = 2393.31 +/- 69.1408 -69.9609 +68.5385 (Minos)
p1 = 6.19428 +/- 0.0759759 -0.0752988 +0.0766601 (Minos)
p2 = 1.68865 +/- 0.0464989 -0.0447029 +0.048734 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2393.31 6.19428 1.68865]
histOutputAgorep14
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 30.8966
NDf = 0
Edm = 2.30646e-08
NCalls = 113
p0 = 2392.11 +/- 71.1208 -72.1391 +70.3889 (Minos)
p1 = 6.18428 +/- 0.0763225 -0.075653 +0.0769982 (Minos)
p2 = 1.691 +/- 0.0484382 -0.0464172 +0.0509987 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2392.11 6.18428 1.691]
histOutputAgorep15
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 28.9404
NDf = 0
Edm = 1.4588e-05
NCalls = 101
p0 = 2389.62 +/- 73.0541 -74.5329 +71.9485 (Minos)
p1 = 6.17486 +/- 0.0766853 -0.0758192 +0.0775703 (Minos)
p2 = 1.69423 +/- 0.0503883 -0.0479924 +0.0534617 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2389.62 6.17486 1.69423]
histOutputAgorep16
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 27.2043
NDf = 0
Edm = 2.00745e-08
NCalls = 108
p0 = 2385.53 +/- 74.9256 -76.3679 +73.9118 (Minos)
p1 = 6.1665 +/- 0.0770866 -0.076446 +0.0777358 (Minos)
p2 = 1.69854 +/- 0.052353 -0.0498695 +0.0555839 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2385.53 6.1665 1.69854]
histOutputAgorep17
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 25.6479
NDf = 0
Edm = 7.49305e-06
NCalls = 98
p0 = 2381.07 +/- 76.6574 -78.4476 +75.3746 (Minos)
p1 = 6.15856 +/- 0.0775109 -0.0769802 +0.0780537 (Minos)
p2 = 1.70313 +/- 0.0542145 -0.0514552 +0.0578389 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2381.07 6.15856 1.70313]
histOutputAgorep18
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 24.2424
NDf = 0
Edm = 3.51473e-08
NCalls = 109
p0 = 2375.84 +/- 78.2406 -79.9679 +77.018 (Minos)
p1 = 6.15095 +/- 0.0779577 -0.0773953 +0.0785384 (Minos)
p2 = 1.70817 +/- 0.055982 -0.0531101 +0.0597242 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2375.84 6.15095 1.70817]
histOutputAgorep19
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 22.9669
NDf = 0
Edm = 7.96352e-08
NCalls = 107
p0 = 2370.62 +/- 79.5865 -81.3434 +78.3158 (Minos)
p1 = 6.14388 +/- 0.0784175 -0.0779004 +0.0789615 (Minos)
p2 = 1.71318 +/- 0.0575408 -0.0545378 +0.061404 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2370.62 6.14388 1.71318]
histOutputAgorep20
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 21.8061
NDf = 0
Edm = 6.90996e-08
NCalls = 109
p0 = 2365.46 +/- 80.6738 -82.3344 +79.4269 (Minos)
p1 = 6.13724 +/- 0.0788756 -0.0783743 +0.0794136 (Minos)
p2 = 1.71809 +/- 0.0588664 -0.0558061 +0.0627092 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2365.46 6.13724 1.71809]
histOutputAgorep21
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 20.7482
NDf = 0
Edm = 3.60576e-08
NCalls = 108
p0 = 2360.59 +/- 81.4757 -82.9796 +80.2862 (Minos)
p1 = 6.13112 +/- 0.0793164 -0.0788597 +0.0798191 (Minos)
p2 = 1.72273 +/- 0.05992 -0.0568573 +0.0636503 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2360.59 6.13112 1.72273]
histOutputAgorep22
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 19.7836
NDf = 0
Edm = 1.67419e-07
NCalls = 106
p0 = 2356.12 +/- 81.9995 -83.2875 +80.9154 (Minos)
p1 = 6.1255 +/- 0.0797259 -0.0793155 +0.0801896 (Minos)
p2 = 1.72699 +/- 0.0606994 -0.0576875 +0.0642428 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2356.12 6.1255 1.72699]
histOutputAgorep23
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 18.9042
NDf = 0
Edm = 7.23017e-07
NCalls = 108
p0 = 2352.05 +/- 82.2803 -83.2445 +81.4023 (Minos)
p1 = 6.12036 +/- 0.0800955 -0.0797097 +0.0805351 (Minos)
p2 = 1.73088 +/- 0.0612334 -0.0583585 +0.0644839 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2352.05 6.12036 1.73088]
histOutputAgorep24
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 18.1027
NDf = 0
Edm = 1.12284e-06
NCalls = 108
p0 = 2348.51 +/- 82.3552 -83.0365 +81.6583 (Minos)
p1 = 6.11569 +/- 0.0804179 -0.08004 +0.0808515 (Minos)
p2 = 1.73431 +/- 0.0615464 -0.0588188 +0.0645108 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2348.51 6.11569 1.73431]
histOutputAgorep25
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 17.3722
NDf = 0
Edm = 2.4526e-07
NCalls = 109
p0 = 2345.48 +/- 82.275 -82.7168 +81.7297 (Minos)
p1 = 6.11149 +/- 0.0806932 -0.0803289 +0.0811141 (Minos)
p2 = 1.73731 +/- 0.0616846 -0.0591194 +0.0643657 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2345.48 6.11149 1.73731]
histOutputAgorep30
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 14.582
NDf = 0
Edm = 1.0438e-07
NCalls = 108
p0 = 2335.86 +/- 80.7919 -80.3191 +81.0186 (Minos)
p1 = 6.09634 +/- 0.0815191 -0.0811746 +0.0819027 (Minos)
p2 = 1.74717 +/- 0.0609954 -0.0592501 +0.0625777 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2335.86 6.09634 1.74717]
histOutputAgorep35
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 12.8045
NDf = 0
Edm = 9.34012e-08
NCalls = 98
p0 = 2331.88 +/- 79.1862 -78.3522 +79.82 (Minos)
p1 = 6.08782 +/- 0.081828 -0.0814668 +0.0822162 (Minos)
p2 = 1.75171 +/- 0.0597919 -0.0585402 +0.0608655 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2331.88 6.08782 1.75171]
histOutputAgorep40
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 11.6423
NDf = 0
Edm = 1.02018e-05
NCalls = 100
p0 = 2330.29 +/- 78.0177 -77.0108 +78.86 (Minos)
p1 = 6.08294 +/- 0.0819366 -0.0815303 +0.0823707 (Minos)
p2 = 1.75394 +/- 0.0588235 -0.057927 +0.0595514 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2330.29 6.08294 1.75394]
histOutputAgorep45
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.8611
NDf = 0
Edm = 1.95297e-08
NCalls = 101
p0 = 2329.76 +/- 77.1875 -76.1446 +78.1275 (Minos)
p1 = 6.08016 +/- 0.0819698 -0.0815648 +0.0823947 (Minos)
p2 = 1.75487 +/- 0.0580944 -0.0572715 +0.0588007 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2329.76 6.08016 1.75487]
histOutputAgorep50
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.3233
NDf = 0
Edm = 4.24138e-08
NCalls = 101
p0 = 2329.68 +/- 76.6148 -75.5421 +77.6124 (Minos)
p1 = 6.07856 +/- 0.0819664 -0.0815546 +0.0823986 (Minos)
p2 = 1.75534 +/- 0.0575737 -0.056851 +0.0582026 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2329.68 6.07856 1.75534]
histOutputAgorep55
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.94615
NDf = 0
Edm = 7.72081e-06
NCalls = 97
p0 = 2329.82 +/- 76.2151 -75.1493 +77.2088 (Minos)
p1 = 6.07766 +/- 0.0819447 -0.0814959 +0.0824191 (Minos)
p2 = 1.75557 +/- 0.0571982 -0.0565833 +0.0577231 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2329.82 6.07766 1.75557]
histOutputAgorep60
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.67782
NDf = 0
Edm = 7.72101e-08
NCalls = 103
p0 = 2329.99 +/- 75.9094 -74.8089 +76.9693 (Minos)
p1 = 6.07728 +/- 0.0819216 -0.0815055 +0.0823564 (Minos)
p2 = 1.75556 +/- 0.0569057 -0.0562971 +0.0574508 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2329.99 6.07728 1.75556]
histOutputAgorep65
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.48518
NDf = 0
Edm = 6.83294e-09
NCalls = 104
p0 = 2330.24 +/- 75.6898 -74.6008 +76.7488 (Minos)
p1 = 6.07715 +/- 0.0818931 -0.0814768 +0.0823276 (Minos)
p2 = 1.75548 +/- 0.0566882 -0.0560956 +0.0572289 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2330.24 6.07715 1.75548]
histOutputAgorep70
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.34625
NDf = 0
Edm = 7.40461e-09
NCalls = 104
p0 = 2330.47 +/- 75.5249 -74.4325 +76.5957 (Minos)
p1 = 6.0772 +/- 0.0818641 -0.0814448 +0.0823013 (Minos)
p2 = 1.75538 +/- 0.0565213 -0.055952 +0.0570466 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2330.47 6.0772 1.75538]
histOutputAgorep75
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.24613
NDf = 0
Edm = 3.60651e-06
NCalls = 90
p0 = 2330.57 +/- 75.3959 -74.1813 +76.6033 (Minos)
p1 = 6.07746 +/- 0.081838 -0.0815121 +0.0821732 (Minos)
p2 = 1.75534 +/- 0.0563932 -0.0559364 +0.0568083 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2330.57 6.07746 1.75534]
histOutputAgorep80
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.17445
NDf = 0
Edm = 2.32805e-08
NCalls = 108
p0 = 2330.9 +/- 75.3033 -74.2008 +76.397 (Minos)
p1 = 6.0776 +/- 0.0818089 -0.0813833 +0.0822519 (Minos)
p2 = 1.7551 +/- 0.0562883 -0.0557547 +0.0567898 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2330.9 6.0776 1.7551]
histOutputAgorep85
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.12391
NDf = 0
Edm = 9.66783e-06
NCalls = 97
p0 = 2330.9 +/- 75.2171 -73.9318 +76.523 (Minos)
p1 = 6.07794 +/- 0.081787 -0.0814075 +0.082173 (Minos)
p2 = 1.75504 +/- 0.0562036 -0.0557696 +0.0566162 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2330.9 6.07794 1.75504]
histOutputAgorep90
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.0892
NDf = 0
Edm = 3.80927e-08
NCalls = 105
p0 = 2331.28 +/- 75.1734 -74.086 +76.2571 (Minos)
p1 = 6.07819 +/- 0.0817594 -0.0813249 +0.0822118 (Minos)
p2 = 1.75479 +/- 0.0561416 -0.0556206 +0.0566372 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2331.28 6.07819 1.75479]
histOutputAgorep95
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.06651
NDf = 0
Edm = 7.61055e-07
NCalls = 104
p0 = 2331.46 +/- 75.1317 -74.0559 +76.2015 (Minos)
p1 = 6.07855 +/- 0.0817372 -0.0813378 +0.0821527 (Minos)
p2 = 1.75465 +/- 0.0560907 -0.0555865 +0.0565685 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2331.46 6.07855 1.75465]
histOutputAgorep100
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.05296
NDf = 0
Edm = 1.02031e-05
NCalls = 98
p0 = 2331.55 +/- 75.0899 -73.9714 +76.2204 (Minos)
p1 = 6.07861 +/- 0.0817142 -0.0810451 +0.0824076 (Minos)
p2 = 1.75449 +/- 0.0560425 -0.0555293 +0.0565429 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2331.55 6.07861 1.75449]
histOutputAgorep150
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.12732
NDf = 0
Edm = 3.08834e-10
NCalls = 109
p0 = 2332.43 +/- 75.0022 -73.894 +76.1246 (Minos)
p1 = 6.08166 +/- 0.0815863 -0.0811589 +0.0820312 (Minos)
p2 = 1.75329 +/- 0.0558779 -0.0554108 +0.0563366 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2332.43 6.08166 1.75329]
histOutputAgorep200
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.27969
NDf = 0
Edm = 1.61096e-07
NCalls = 111
p0 = 2332.8 +/- 75.0377 -73.8959 +76.1994 (Minos)
p1 = 6.08337 +/- 0.0815398 -0.0811157 +0.0819812 (Minos)
p2 = 1.75248 +/- 0.0558583 -0.055423 +0.0562881 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2332.8 6.08337 1.75248]
histOutputAgorep250
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.41695
NDf = 0
Edm = 2.24257e-08
NCalls = 110
p0 = 2333.07 +/- 75.0859 -73.9573 +76.2347 (Minos)
p1 = 6.08433 +/- 0.0815362 -0.081113 +0.081977 (Minos)
p2 = 1.75185 +/- 0.0558629 -0.055419 +0.0563031 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.07 6.08433 1.75185]
histOutputAgorep300
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.53327
NDf = 0
Edm = 6.87797e-08
NCalls = 110
p0 = 2333.25 +/- 75.1273 -73.9901 +76.2857 (Minos)
p1 = 6.08484 +/- 0.081555 -0.0811333 +0.0819944 (Minos)
p2 = 1.75138 +/- 0.0558714 -0.0554358 +0.0563036 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.25 6.08484 1.75138]
histOutputAgorep350
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.63207
NDf = 0
Edm = 1.47821e-07
NCalls = 110
p0 = 2333.38 +/- 75.1608 -74.0152 +76.328 (Minos)
p1 = 6.08509 +/- 0.0815844 -0.0811639 +0.0820226 (Minos)
p2 = 1.75103 +/- 0.0558793 -0.0554515 +0.0563036 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.38 6.08509 1.75103]
histOutputAgorep400
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.71683
NDf = 0
Edm = 2.36995e-07
NCalls = 110
p0 = 2333.48 +/- 75.1872 -74.0346 +76.3617 (Minos)
p1 = 6.08518 +/- 0.0816173 -0.0811977 +0.0820549 (Minos)
p2 = 1.75075 +/- 0.0558861 -0.0554647 +0.0563037 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.48 6.08518 1.75075]
histOutputAgorep450
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.7902
NDf = 0
Edm = 3.13329e-07
NCalls = 110
p0 = 2333.55 +/- 75.2079 -74.0499 +76.388 (Minos)
p1 = 6.08518 +/- 0.0816499 -0.0812307 +0.0820871 (Minos)
p2 = 1.75054 +/- 0.0558917 -0.0554751 +0.0563043 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.55 6.08518 1.75054]
histOutputAgorep500
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.85409
NDf = 0
Edm = 3.63948e-07
NCalls = 110
p0 = 2333.59 +/- 75.2241 -74.0622 +76.4083 (Minos)
p1 = 6.08512 +/- 0.08168 -0.081261 +0.0821172 (Minos)
p2 = 1.75037 +/- 0.0558964 -0.055483 +0.0563054 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.59 6.08512 1.75037]
histOutputAgorep550
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.90995
NDf = 0
Edm = 3.86269e-07
NCalls = 110
p0 = 2333.62 +/- 75.2368 -74.0722 +76.4239 (Minos)
p1 = 6.08502 +/- 0.0817068 -0.0812877 +0.082144 (Minos)
p2 = 1.75025 +/- 0.0559002 -0.0554888 +0.0563071 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.62 6.08502 1.75025]
histOutputAgorep600
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 9.95893
NDf = 0
Edm = 3.84226e-07
NCalls = 110
p0 = 2333.62 +/- 75.247 -74.0806 +76.4358 (Minos)
p1 = 6.08492 +/- 0.08173 -0.0813106 +0.0821674 (Minos)
p2 = 1.75015 +/- 0.0559032 -0.0554929 +0.056309 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.62 6.08492 1.75015]
histOutputAgorep650
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.0019
NDf = 0
Edm = 3.64514e-07
NCalls = 110
p0 = 2333.62 +/- 75.2551 -74.0877 +76.445 (Minos)
p1 = 6.08481 +/- 0.0817498 -0.0813301 +0.0821875 (Minos)
p2 = 1.75008 +/- 0.0559055 -0.0554956 +0.0563109 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.62 6.08481 1.75008]
histOutputAgorep700
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.0397
NDf = 0
Edm = 3.34042e-07
NCalls = 110
p0 = 2333.61 +/- 75.2616 -74.0938 +76.4522 (Minos)
p1 = 6.0847 +/- 0.0817665 -0.0813465 +0.0822047 (Minos)
p2 = 1.75002 +/- 0.0559072 -0.0554971 +0.0563128 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.61 6.0847 1.75002]
histOutputAgorep750
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.0729
NDf = 0
Edm = 2.9859e-07
NCalls = 110
p0 = 2333.59 +/- 75.2671 -74.0991 +76.4578 (Minos)
p1 = 6.0846 +/- 0.0817807 -0.0813602 +0.0822193 (Minos)
p2 = 1.74997 +/- 0.0559083 -0.0554977 +0.0563144 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.59 6.0846 1.74997]
histOutputAgorep800
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.1021
NDf = 0
Edm = 2.62333e-07
NCalls = 110
p0 = 2333.57 +/- 75.2717 -74.1039 +76.4624 (Minos)
p1 = 6.0845 +/- 0.0817927 -0.0813717 +0.0822317 (Minos)
p2 = 1.74993 +/- 0.0559088 -0.0554975 +0.0563156 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.57 6.0845 1.74993]
histOutputAgorep850
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.1277
NDf = 0
Edm = 2.27907e-07
NCalls = 110
p0 = 2333.54 +/- 75.2758 -74.1082 +76.4661 (Minos)
p1 = 6.08441 +/- 0.0818029 -0.0813815 +0.0822424 (Minos)
p2 = 1.74989 +/- 0.0559088 -0.0554968 +0.0563164 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.54 6.08441 1.74989]
histOutputAgorep900
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.1502
NDf = 0
Edm = 1.96724e-07
NCalls = 110
p0 = 2333.52 +/- 75.2793 -74.1122 +76.4694 (Minos)
p1 = 6.08433 +/- 0.0818116 -0.0813897 +0.0822515 (Minos)
p2 = 1.74986 +/- 0.0559084 -0.0554956 +0.0563168 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.52 6.08433 1.74986]
histOutputAgorep950
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.17
NDf = 0
Edm = 1.69332e-07
NCalls = 110
p0 = 2333.5 +/- 75.2825 -74.1158 +76.4722 (Minos)
p1 = 6.08426 +/- 0.081819 -0.0813967 +0.0822593 (Minos)
p2 = 1.74983 +/- 0.0559076 -0.0554941 +0.0563169 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.5 6.08426 1.74983]
histOutputAgorep1000
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 10.1874
NDf = 0
Edm = 1.45744e-07
NCalls = 110
p0 = 2333.49 +/- 75.2855 -74.1192 +76.4748 (Minos)
p1 = 6.08419 +/- 0.0818254 -0.0814028 +0.0822661 (Minos)
p2 = 1.7498 +/- 0.0559066 -0.0554924 +0.0563165 (Minos)
fcn flag=0: npar=3 gin=0 par=[ 2333.49 6.08419 1.7498]
"inversion",0.369615,16,0.989,0.687185,13,0.777954,0.843936,6.07588,"±",0.0843195,1.76119,"±",0.0622688
"template",0.220299,16,0.99951,0.640838,13,0.821395,0.834923,6.0667,"±",0.0821746,1.77123,"±",0.0599625
"template+area",0.237251,16,0.999212,0.643779,13,0.818749,0.834905,6.06894,"±",0.0818049,1.77403,"±",0.0597178
"Tikhonov+area",0.567834,16,0.909857,1.70353,13,0.0531425,0.666522,6.01392,"±",0.0830555,1.88219,"±",0.0523609
"min(rhomax)",1.30367,16,0.184004,2.80306,13,0.000506874,0.656519,5.94814,"±",0.0813071,1.96606,"±",0.0470327
"bin-by-bin",4.925,16,2.7358e-10,1.46604,13,0.121301,0,5.90214,"±",0.0816763,2.07111,"±",0.0460756
"iterative, N=0",319.06,16,0,15.2271,13,3.5961e-35,0.849481,5.62399,"±",0.0532142,2.24364,"±",0.020155
"iterative, N=1",111.402,16,0,9.09211,13,4.5692e-19,0.817166,5.8135,"±",0.0699476,2.27527,"±",0.0262139
"iterative, N=2",62.992,16,2.31021e-204,7.97086,13,3.28626e-16,0.785917,5.89902,"±",0.0771589,2.27417,"±",0.0289737
"iterative, N=3",41.9938,16,1.22802e-132,7.44678,13,6.87329e-15,0.757593,5.93729,"±",0.0805321,2.25889,"±",0.0305052
"iterative, N=4",30.3819,16,2.83616e-93,7.04813,13,6.82474e-14,0.733895,5.9568,"±",0.0823367,2.2385,"±",0.0316079
"iterative, N=5",23.0971,16,8.57418e-69,6.68416,13,5.46969e-13,0.715893,5.96958,"±",0.0835144,2.2164,"±",0.0326294
"iterative, N=6",18.1625,16,2.24797e-52,6.33352,13,4.0031e-12,0.703607,5.98002,"±",0.0844752,2.1935,"±",0.0337655
"iterative, N=7",14.6425,16,8.54613e-41,5.99028,13,2.76604e-11,0.696048,5.99052,"±",0.0854688,2.16994,"±",0.0352059
"iterative, N=8",12.0361,16,2.49723e-32,5.65178,13,1.83018e-10,0.691821,6.0021,"±",0.0867252,2.14524,"±",0.0372357
"iterative, N=9",10.0504,16,5.68598e-26,5.31529,13,1.17568e-09,0.6897,6.01603,"±",0.0887007,2.11836,"±",0.0405096
"iterative, N=10",8.50328,16,4.25835e-21,3.27029,13,5.39585e-05,0.688825,6.23206,"±",0.0749671,1.6911,"±",0.0410841
"iterative, N=11",7.27553,16,2.68544e-17,2.98421,13,0.000215602,0.688651,6.21773,"±",0.0753106,1.6883,"±",0.0428048
"iterative, N=12",6.28618,16,2.69801e-14,2.74799,13,0.000654782,0.688847,6.20537,"±",0.0756438,1.68766,"±",0.0446165
"iterative, N=13",5.47848,16,6.74681e-12,2.5485,13,0.00162879,0.689213,6.19428,"±",0.0759759,1.68865,"±",0.0464989
"iterative, N=14",4.81161,16,5.78501e-10,2.37666,13,0.00349095,0.68963,6.18428,"±",0.0763225,1.691,"±",0.0484382
"iterative, N=15",4.25556,16,2.15151e-08,2.22618,13,0.00667463,0.690033,6.17486,"±",0.0766853,1.69423,"±",0.0503883
"iterative, N=16",3.78782,16,4.13727e-07,2.09264,13,0.0116635,0.690388,6.1665,"±",0.0770866,1.69854,"±",0.052353
"iterative, N=17",3.3913,16,4.70074e-06,1.97292,13,0.0189482,0.690682,6.15856,"±",0.0775109,1.70313,"±",0.0542145
"iterative, N=18",3.05275,16,3.49724e-05,1.8648,13,0.0289726,0.690913,6.15095,"±",0.0779577,1.70817,"±",0.055982
"iterative, N=19",2.76185,16,0.000184579,1.76668,13,0.0420756,0.691087,6.14388,"±",0.0784175,1.71318,"±",0.0575408
"iterative, N=20",2.51042,16,0.000736128,1.67739,13,0.0584402,0.691212,6.13724,"±",0.0788756,1.71809,"±",0.0588664
"iterative, N=21",2.29192,16,0.00233267,1.59601,13,0.0780608,0.691301,6.13112,"±",0.0793164,1.72273,"±",0.05992
"iterative, N=22",2.10111,16,0.00611486,1.52182,13,0.100738,0.691362,6.1255,"±",0.0797259,1.72699,"±",0.0606994
"iterative, N=23",1.93369,16,0.0136996,1.45417,13,0.126103,0.691407,6.12036,"±",0.0800955,1.73088,"±",0.0612334
"iterative, N=24",1.78617,16,0.0269367,1.39251,13,0.15366,0.691446,6.11569,"±",0.0804179,1.73431,"±",0.0615464
"iterative, N=25",1.65566,16,0.0475044,1.33632,13,0.182842,0.691486,6.11149,"±",0.0806932,1.73731,"±",0.0616846
"iterative, N=30",1.18689,16,0.269169,1.12169,13,0.334167,0.691931,6.09634,"±",0.0815191,1.74717,"±",0.0609954
"iterative, N=35",0.908322,16,0.559059,0.984963,13,0.463022,0.69311,6.08782,"±",0.081828,1.75171,"±",0.0597919
"iterative, N=40",0.732607,16,0.762905,0.895559,13,0.557174,0.695148,6.08294,"±",0.0819366,1.75394,"±",0.0588235
"iterative, N=45",0.616619,16,0.873537,0.835466,13,0.622456,0.69795,6.08016,"±",0.0819698,1.75487,"±",0.0580944
"iterative, N=50",0.537326,16,0.92907,0.794103,13,0.667304,0.701338,6.07856,"±",0.0819664,1.75534,"±",0.0575737
"iterative, N=55",0.481626,16,0.95714,0.765088,13,0.698332,0.705124,6.07766,"±",0.0819447,1.75557,"±",0.0571982
"iterative, N=60",0.441666,16,0.971955,0.744448,13,0.720044,0.709144,6.07728,"±",0.0819216,1.75556,"±",0.0569057
"iterative, N=65",0.412534,16,0.980214,0.729629,13,0.735393,0.713267,6.07715,"±",0.0818931,1.75548,"±",0.0566882
"iterative, N=70",0.391042,16,0.985078,0.718942,13,0.746319,0.717396,6.0772,"±",0.0818641,1.75538,"±",0.0565213
"iterative, N=75",0.375056,16,0.988086,0.71124,13,0.754111,0.72146,6.07746,"±",0.081838,1.75534,"±",0.0563932
"iterative, N=80",0.363112,16,0.990023,0.705727,13,0.759644,0.725414,6.0776,"±",0.0818089,1.7551,"±",0.0562883
"iterative, N=85",0.354178,16,0.991311,0.701839,13,0.763522,0.729226,6.07794,"±",0.081787,1.75504,"±",0.0562036
"iterative, N=90",0.347516,16,0.992188,0.699169,13,0.766174,0.732878,6.07819,"±",0.0817594,1.75479,"±",0.0561416
"iterative, N=95",0.342588,16,0.992792,0.697424,13,0.767902,0.736361,6.07855,"±",0.0817372,1.75465,"±",0.0560907
"iterative, N=100",0.338996,16,0.99321,0.696382,13,0.768932,0.739673,6.07861,"±",0.0817142,1.75449,"±",0.0560425
"iterative, N=150",0.335268,16,0.993624,0.702102,13,0.763261,0.764528,6.08166,"±",0.0815863,1.75329,"±",0.0558779
"iterative, N=200",0.346127,16,0.992362,0.713822,13,0.751507,0.779298,6.08337,"±",0.0815398,1.75248,"±",0.0558583
"iterative, N=250",0.356914,16,0.990931,0.724381,13,0.740775,0.788781,6.08433,"±",0.0815362,1.75185,"±",0.0558629
"iterative, N=300",0.36599,16,0.989579,0.733328,13,0.731582,0.795305,6.08484,"±",0.081555,1.75138,"±",0.0558714
"iterative, N=350",0.373445,16,0.988362,0.740928,13,0.723709,0.800037,6.08509,"±",0.0815844,1.75103,"±",0.0558793
"iterative, N=400",0.379583,16,0.987284,0.747449,13,0.71691,0.803611,6.08518,"±",0.0816173,1.75075,"±",0.0558861
"iterative, N=450",0.384672,16,0.986336,0.753092,13,0.710994,0.806396,6.08518,"±",0.0816499,1.75054,"±",0.0558917
"iterative, N=500",0.388923,16,0.985505,0.758007,13,0.705821,0.808621,6.08512,"±",0.08168,1.75037,"±",0.0558964
"iterative, N=550",0.3925,16,0.984779,0.762304,13,0.701281,0.810434,6.08502,"±",0.0817068,1.75025,"±",0.0559002
"iterative, N=600",0.395528,16,0.984143,0.766072,13,0.697289,0.811937,6.08492,"±",0.08173,1.75015,"±",0.0559032
"iterative, N=650",0.398104,16,0.983588,0.76938,13,0.693776,0.813199,6.08481,"±",0.0817498,1.75008,"±",0.0559055
"iterative, N=700",0.400304,16,0.983103,0.772286,13,0.690683,0.814272,6.0847,"±",0.0817665,1.75002,"±",0.0559072
"iterative, N=750",0.40219,16,0.982678,0.77484,13,0.687961,0.815194,6.0846,"±",0.0817807,1.74997,"±",0.0559083
"iterative, N=800",0.403812,16,0.982307,0.777084,13,0.685565,0.815995,6.0845,"±",0.0817927,1.74993,"±",0.0559088
"iterative, N=850",0.40521,16,0.981983,0.779056,13,0.683458,0.816697,6.08441,"±",0.0818029,1.74989,"±",0.0559088
"iterative, N=900",0.406417,16,0.9817,0.780788,13,0.681605,0.817318,6.08433,"±",0.0818116,1.74986,"±",0.0559084
"iterative, N=950",0.407463,16,0.981452,0.78231,13,0.679975,0.817871,6.08426,"±",0.081819,1.74983,"±",0.0559076
"iterative, N=1000",0.408369,16,0.981235,0.783648,13,0.678541,0.818368,6.08419,"±",0.0818254,1.7498,"±",0.0559066
#include <iostream>
#include <cmath>
#include <map>
#include <TMath.h>
#include <TCanvas.h>
#include <TStyle.h>
#include <TGraph.h>
#include <TGraphErrors.h>
#include <TFile.h>
#include <TROOT.h>
#include <TText.h>
#include <TLine.h>
#include <TLegend.h>
#include <TH1.h>
#include <TF1.h>
#include <TFitter.h>
#include <TMatrixD.h>
#include <TMatrixDSym.h>
#include <TVectorD.h>
#include <TFitResult.h>
#include <TRandom3.h>
#include "TUnfoldDensity.h"
using std::vector, std::pair, std::cout;
// #define PRINT_MATRIX_L
#define TEST_INPUT_COVARIANCE
void CreateHistogramCopies(TH1 *h[3],TUnfoldBinning const *binning);
TH2 *AddOverflowXY(TH2 *h,double widthX,double widthY);
void DrawOverflowX(TH1 *h,double posy);
void DrawOverflowY(TH1 *h,double posx);
double const kLegendFontSize=0.05;
int kNbinC=0;
char const *text=nullptr,double tau=0.0,vector<double> const *r=nullptr,
TF1 *f=nullptr);
vector<pair<TF1*,vector<double> > > const *table);
TFitResultPtr DoFit(TH1 *h,TH2 *rho,TH1 *truth,char const *text,
vector<pair<TF1*,vector<double> > > &table,int niter=0);
void GetNiterGraphs(int iFirst,int iLast,vector<pair<TF1*,
vector<double> > > const &table,int color,
TGraph *graph[4],int style);
void GetNiterHist(int ifit,vector<pair<TF1*,vector<double> > > const &table,
TH1 *hist[4],int color,int style,int fillStyle);
#ifdef WITH_IDS
#endif
{
// switch on histogram errors
g_rnd=new TRandom3(4711);
//==============================================
// step 1 : open output file
TFile *outputFile=new TFile("testUnfold7_results.root","recreate");
//==============================================
// step 2 : read binning schemes and input histograms
TFile *inputFile=new TFile("testUnfold7_histograms.root");
outputFile->cd();
inputFile->GetObject("fine",fineBinning);
inputFile->GetObject("coarse",coarseBinning);
cout<<"problem to read binning schemes\n";
}
// save binning schemes to output file
fineBinning->Write();
coarseBinning->Write();
// read histograms
#define READ(TYPE,binning,name) \
TYPE *name[3]; inputFile->GetObject(#name,name[0]); \
name[0]->Write(); \
if(!name[0]) cout<<"Error reading " #name "\n"; \
CreateHistogramCopies(name,binning);
outputFile->cd();
//TH2 *histRhoCLCurve;
double tauMin=1.e-4;
double tauMax=1.e-1;
double fBgr=1.0; // 0.2/0.25;
double biasScale=1.0;
//double tauC;
{
tunfoldC->SubtractBackground(histMcbgrRecC[0],"BGR",fBgr,0.0);
tunfoldC->DoUnfold(0.);
histOutputCtau0[0]=tunfoldC->GetOutput("histOutputCtau0");
histRhoCtau0=tunfoldC->GetRhoIJtotal("histRhoCtau0");
tunfoldC->ScanLcurve(50,tauMin,tauMax,nullptr);
/* tauC= */tunfoldC->GetTau();
//tunfoldC->ScanTau(50,1.E-7,1.E-1,0,TUnfoldDensity::kEScanTauRhoAvg);
histOutputCLCurve[0]=tunfoldC->GetOutput("histOutputCLCurve");
/* histRhoCLCurve= */tunfoldC->GetRhoIJtotal("histRhoCLCurve");
histProbC=tunfoldC->GetProbabilityMatrix("histProbC",";P_T(gen);P_T(rec)");
}
//TH2 *histRhoFLCurve;
tauMin=3.E-4;
tauMax=3.E-2;
//double tauF;
{
tunfoldF->SubtractBackground(histMcbgrRecF[0],"BGR",fBgr,0.0);
tunfoldF->DoUnfold(0.);
histOutputFtau0[0]=tunfoldF->GetOutput("histOutputFtau0");
histRhoFtau0=tunfoldF->GetRhoIJtotal("histRhoFtau0");
tunfoldF->ScanLcurve(50,tauMin,tauMax,nullptr);
//tunfoldF->DoUnfold(tauC);
/* tauF= */tunfoldF->GetTau();
//tunfoldF->ScanTau(50,1.E-7,1.E-1,0,TUnfoldDensity::kEScanTauRhoAvg);
histOutputFLCurve[0]=tunfoldF->GetOutput("histOutputFLCurve");
/* histRhoFLCurve= */tunfoldF->GetRhoIJtotal("histRhoFLCurve");
histProbF=tunfoldF->GetProbabilityMatrix("histProbF",";P_T(gen);P_T(rec)");
}
TSpline *rhoScan=nullptr;
double tauFA,tauFArho;
{
tunfoldFA->SubtractBackground(histMcbgrRecF[0],"BGR",fBgr,0.0);
tunfoldFA->DoUnfold(0.);
histOutputFAtau0[0]=tunfoldFA->GetOutput("histOutputFAtau0");
histRhoFAtau0=tunfoldFA->GetRhoIJtotal("histRhoFAtau0");
tauFArho=tunfoldFA->GetTau();
histOutputFArho[0]=tunfoldFA->GetOutput("histOutputFArho");
histRhoFArho=tunfoldFA->GetRhoIJtotal("histRhoFArho");
tauFA=tunfoldFA->GetTau();
histOutputFALCurve[0]=tunfoldFA->GetOutput("histOutputFALCurve");
histRhoFALCurve=tunfoldFA->GetRhoIJtotal("histRhoFALCurve");
}
lCurve->Write();
logTauX->Write();
logTauY->Write();
double widthC=coarseBinning->GetBinSize(histProbC->GetNbinsY()+1);
double widthF=fineBinning->GetBinSize(histProbF->GetNbinsY()+1);
// efficiency
TH1 *histEfficiencyC=histProbCO->ProjectionX("histEfficiencyC");
kNbinC=histProbCO->GetNbinsX();
// reconstructed quantities with overflow (coarse binning)
// MC: add signal and bgr
histMcbgrRecCO->Scale(fBgr);
TH1 *histMcRecCO=(TH1 *)histMcsigRecCO->Clone("histMcRecC0");
//TH1 *histDataRecCNO=AddOverflowX(histDataRecC[1],widthC);
histMcbgrRecFO->Scale(fBgr);
TH1 *histMcRecFO=(TH1 *)histMcsigRecFO->Clone("histMcRecF0");
// truth level with overflow
// unfolding result with overflow
//TH1 *histOutputCLCurveO=AddOverflowX(histOutputCLCurve[2],widthC);
//TH2 *histRhoCLCurveO=AddOverflowXY(histRhoCLCurve,widthC,widthC);
//TH1 *histOutputFLCurveO=AddOverflowX(histOutputFLCurve[2],widthC);
//TH2 *histRhoFLCurveO=AddOverflowXY(histRhoFLCurve,widthC,widthC);
// bin-by-bin
TH2 *histRhoBBBO=(TH2 *)histRhoCtau0O->Clone("histRhoBBBO");
for(int i=1;i<=histRhoBBBO->GetNbinsX();i++) {
for(int j=1;j<=histRhoBBBO->GetNbinsX();j++) {
histRhoBBBO->SetBinContent(i,j,(i==j)?1.:0.);
}
}
TH1 *histDataBgrsub=(TH1 *)histDataRecCO->Clone("histDataBgrsub");
TH1 *histOutputBBBO=(TH1 *)histDataBgrsub->Clone("histOutputBBBO");
// iterative
int niter=1000;
cout<<"maximum number of iterations: "<<niter<<"\n";
histOutputAgo.push_back((TH1*)histMcsigGenO->Clone("histOutputAgo-1"));
histOutputAgorep.push_back((TH1*)histMcsigGenO->Clone("histOutputAgorep-1"));
histRhoAgo.push_back((TH2*)histRhoBBBO->Clone("histRhoAgo-1"));
histRhoAgorep.push_back((TH2*)histRhoBBBO->Clone("histRhoAgorep-1"));
nIter.push_back(-1);
int nx=histProbCO->GetNbinsX();
int ny=histProbCO->GetNbinsY();
for(int i=0;i<nx;i++) {
double epsilonI=0.;
for(int j=0;j<ny;j++) {
epsilonI+= histProbCO->GetBinContent(i+1,j+1);
}
for(int j=0;j<ny;j++) {
double aji=histProbCO->GetBinContent(i+1,j+1);
A(j,i)=aji;
}
}
for(int i=0;i<nx;i++) {
(histOutputAgo[0]->GetBinError(i+1)
*histOutputAgo[0]->GetXaxis()->GetBinWidth(i+1),2.);
}
for(int i=0;i<ny;i++) {
(histDataRecCO->GetBinError(i+1)
*histDataRecCO->GetXaxis()->GetBinWidth(i+1),2.);
}
#define NREPLICA 300
for(int nr=0;nr<NREPLICA;nr++) {
x[nr]=new TVectorD(nx);
y[nr]=new TVectorD(ny);
yMb[nr]=new TVectorD(ny);
yErr[nr]=new TVectorD(ny);
}
for(int i=0;i<nx;i++) {
(*x[0])(i)=histOutputAgo[0]->GetBinContent(i+1)
*histOutputAgo[0]->GetXaxis()->GetBinWidth(i+1);
for(int nr=1;nr<NREPLICA;nr++) {
(*x[nr])(i)=(*x[0])(i);
}
}
for(int i=0;i<ny;i++) {
(*y[0])(i)=histDataRecCO->GetBinContent(i+1)
*histDataRecCO->GetXaxis()->GetBinWidth(i+1);
for(int nr=1;nr<NREPLICA;nr++) {
(*y[nr])(i)=g_rnd->Poisson((*y[0])(i));
(*yErr[nr])(i)=TMath::Sqrt((*y[nr])(i));
}
b(i)=histMcbgrRecCO->GetBinContent(i+1)*
histMcbgrRecCO->GetXaxis()->GetBinWidth(i+1);
for(int nr=0;nr<NREPLICA;nr++) {
(*yMb[nr])(i)=(*y[nr])(i)-b(i);
}
}
for(int iter=0;iter<=niter;iter++) {
if(!(iter %100)) cout<<iter<<"\n";
for(int nr=0;nr<NREPLICA;nr++) {
TVectorD yrec=A*(*x[nr])+b;
for(int j=0;j<ny;j++) {
yOverYrec(j)=(*y[nr])(j)/yrec(j);
}
for(int i=0;i<nx;i++) {
xx(i) = (*x[nr])(i) * f(i);
}
if(nr==0) {
for(int i=0;i<nx;i++) {
for(int j=0;j<ny;j++) {
xdf_dr(i,j) *= (*x[nr])(i);
}
}
for(int j=0;j<ny;j++) {
dr_dxdy(j,nx+j)=1.0/yrec(j);
for(int i=0;i<nx;i++) {
dr_dxdy(j,i)= -yOverYrec(j)/yrec(j)*A(j,i);
}
}
dxy_dxy.SetSub(0,0,xdf_dr*dr_dxdy);
for(int i=0;i<nx;i++) {
dxy_dxy(i,i) +=f(i);
}
for(int i=0;i<ny;i++) {
dxy_dxy(nx+i,nx+i) +=1.0;
}
}
(*x[nr])=xx;
}
if((iter<=25)||
((iter<=100)&&(iter %5==0))||
((iter<=1000)&&(iter %50==0))||
(iter %1000==0)) {
nIter.push_back(iter);
TH1 * h=(TH1*)histOutputAgo[0]->Clone
(TString::Format("histOutputAgo%d",iter));
histOutputAgo.push_back(h);
for(int i=0;i<nx;i++) {
double bw=h->GetXaxis()->GetBinWidth(i+1);
h->SetBinContent(i+1,(*x[0])(i)/bw);
h->SetBinError(i+1,TMath::Sqrt(covAgo(i,i))/bw);
}
TH2 *h2=(TH2*)histRhoAgo[0]->Clone
(TString::Format("histRhoAgo%d",iter));
histRhoAgo.push_back(h2);
for(int i=0;i<nx;i++) {
for(int j=0;j<nx;j++) {
double rho= covAgo(i,j)/TMath::Sqrt(covAgo(i,i)*covAgo(j,j));
if((i!=j)&&(TMath::Abs(rho)>=1.0)) {
cout<<"bad error matrix: iter="<<iter<<"\n";
exit(0);
}
h2->SetBinContent(i+1,j+1,rho);
}
}
// error and correlations from replica analysis
h=(TH1*)histOutputAgo[0]->Clone
(TString::Format("histOutputAgorep%d",iter));
h2=(TH2*)histRhoAgo[0]->Clone
(TString::Format("histRhoAgorep%d",iter));
histOutputAgorep.push_back(h);
histRhoAgorep.push_back(h2);
TVectorD mean(nx);
double w=1./(NREPLICA-1.);
for(int nr=1;nr<NREPLICA;nr++) {
mean += w* *x[nr];
}
for(int nr=1;nr<NREPLICA;nr++) {
//TMatrixD dx= (*x)-mean;
for(int i=0;i<nx;i++) {
dx(i,0)= (*x[nr])(i)-(*x[0])(i);
}
}
for(int i=0;i<nx;i++) {
double bw=h->GetXaxis()->GetBinWidth(i+1);
h->SetBinContent(i+1,(*x[0])(i)/bw);
h->SetBinError(i+1,TMath::Sqrt(covAgorep(i,i))/bw);
// cout<<i<<" "<<(*x[0])(i)/bw<<" +/-"<<TMath::Sqrt(covAgorep(i,i))/bw<<" "<<TMath::Sqrt(covAgo(i,i))/bw<<"\n";
}
for(int i=0;i<nx;i++) {
for(int j=0;j<nx;j++) {
double rho= covAgorep(i,j)/
if((i!=j)&&(TMath::Abs(rho)>=1.0)) {
cout<<"bad error matrix: iter="<<iter<<"\n";
exit(0);
}
h2->SetBinContent(i+1,j+1,rho);
}
}
}
}
#ifdef WITH_IDS
// IDS Malaescu
int niterIDS=100;
cout<<"IDS number of iterations: "<<niterIDS<<"\n";
for(int nr=0;nr<NREPLICA;nr++) {
}
for(int iy=0;iy<ny;iy++) {
for(int ix=0;ix<nx;ix++) {
A_IDS(iy,ix)=histMcsigGenRecC[0]->GetBinContent(ix+1,iy+1);
}
}
double lambdaL=0.;
Double_t lambdaUmin = 1.0000002;
Double_t lambdaMmin = 0.0000001;
Double_t lambdaS = 0.000001;
histOutputIDS.push_back((TH1*)histOutputAgo[0]->Clone("histOutputIDS-1"));
histRhoIDS.push_back((TH2*)histRhoAgo[0]->Clone("histRhoIDS-1"));
histOutputIDS.push_back((TH1*)histOutputAgo[0]->Clone("histOutputIDS0"));
histRhoIDS.push_back((TH2*)histRhoAgo[0]->Clone("histRhoIDS0"));
for(int iter=1;iter<=niterIDS;iter++) {
if(!(iter %10)) cout<<iter<<"\n";
for(int nr=0;nr<NREPLICA;nr++) {
if(iter==1) {
} else {
}
}
unsigned ix;
for(ix=0;ix<nIter.size();ix++) {
if(nIter[ix]==iter) break;
}
if(ix<nIter.size()) {
TH1 * h=(TH1*)histOutputIDS[0]->Clone
(TString::Format("histOutputIDS%d",iter));
TH2 *h2=(TH2*)histRhoIDS[0]->Clone
(TString::Format("histRhoIDS%d",iter));
histOutputIDS.push_back(h);
histRhoIDS.push_back(h2);
TVectorD mean(nx);
double w=1./(NREPLICA-1.);
for(int nr=1;nr<NREPLICA;nr++) {
mean += w* (*unfresIDS[nr]);
}
for(int nr=1;nr<NREPLICA;nr++) {
//TMatrixD dx= (*x)-mean;
for(int i=0;i<nx;i++) {
dx(i,0)= (*unfresIDS[nr])(i)-(*unfresIDS[0])(i);
}
}
for(int i=0;i<nx;i++) {
double bw=h->GetXaxis()->GetBinWidth(i+1);
h->SetBinContent(i+1,(*unfresIDS[0])(i)/bw/
histEfficiencyC->GetBinContent(i+1));
h->SetBinError(i+1,TMath::Sqrt(covIDSrep(i,i))/bw/
histEfficiencyC->GetBinContent(i+1));
// cout<<i<<" "<<(*x[0])(i)/bw<<" +/-"<<TMath::Sqrt(covAgorep(i,i))/bw<<" "<<TMath::Sqrt(covAgo(i,i))/bw<<"\n";
}
for(int i=0;i<nx;i++) {
for(int j=0;j<nx;j++) {
double rho= covIDSrep(i,j)/
if((i!=j)&&(TMath::Abs(rho)>=1.0)) {
cout<<"bad error matrix: iter="<<iter<<"\n";
exit(0);
}
h2->SetBinContent(i+1,j+1,rho);
}
}
}
}
#endif
//double NEdSmc=histDataBgrsub->GetSumOfWeights();
TCanvas *c1=new TCanvas("c1","",600,600);
TCanvas *c2sq=new TCanvas("c2sq","",600,600);
c2sq->Divide(1,2);
TCanvas *c2w=new TCanvas("c2w","",600,300);
c2w->Divide(2,1);
TCanvas *c4=new TCanvas("c4","",600,600);
c4->Divide(2,2);
//TCanvas *c3n=new TCanvas("c3n","",600,600);
TPad *subn[3];
//gROOT->SetStyle("xTimes2");
subn[0]= new TPad("subn0","",0.,0.5,1.,1.);
//gROOT->SetStyle("square");
subn[1]= new TPad("subn1","",0.,0.,0.5,0.5);
subn[2]= new TPad("subn2","",0.5,0.0,1.,0.5);
for(int i=0;i<3;i++) {
subn[i]->SetFillStyle(0);
subn[i]->Draw();
}
TCanvas *c3c=new TCanvas("c3c","",600,600);
TPad *subc[3];
//gROOT->SetStyle("xTimes2");
subc[0]= new TPad("sub0","",0.,0.5,1.,1.);
//gROOT->SetStyle("squareCOLZ");
subc[1]= new TPad("sub1","",0.,0.,0.5,0.5);
//gROOT->SetStyle("square");
subc[2]= new TPad("sub2","",0.5,0.0,1.,0.5);
for(int i=0;i<3;i++) {
subc[i]->SetFillStyle(0);
subc[i]->Draw();
}
//=========================== example ==================================
c2w->cd(1);
c2w->cd(2);
c2w->SaveAs("exampleTR.eps");
//=========================== example ==================================
c2w->cd(1);
c2w->cd(2);
c2w->SaveAs("exampleAE.eps");
int iFitInversion=table.size();
DoFit(histOutputCtau0O,histRhoCtau0O,histDataGenO,"inversion",table);
//=========================== inversion ==================================
subc[0]->cd();
&table[table.size()-1].second);
subc[1]->cd();
subc[2]->cd();
c3c->SaveAs("inversion.eps");
DoFit(histOutputFtau0O,histRhoFtau0O,histDataGenO,"template",table);
//=========================== template ==================================
subc[0]->cd();
&table[table.size()-1].second);
subc[1]->cd();
subc[2]->cd();
c3c->SaveAs("template.eps");
DoFit(histOutputFAtau0O,histRhoFAtau0O,histDataGenO,"template+area",table);
//=========================== template+area ==================================
subc[0]->cd();
&table[table.size()-1].second);
subc[1]->cd();
subc[2]->cd();
c3c->SaveAs("templateA.eps");
int iFitFALCurve=table.size();
DoFit(histOutputFALCurveO,histRhoFALCurveO,histDataGenO,"Tikhonov+area",table);
//=========================== template+area+tikhonov =====================
subc[0]->cd();
&table[table.size()-1].second);
subc[1]->cd();
subc[2]->cd();
c3c->SaveAs("lcurveFA.eps");
int iFitFArho=table.size();
DoFit(histOutputFArhoO,histRhoFArhoO,histDataGenO,"min(rhomax)",table);
//=========================== template+area+tikhonov =====================
subc[0]->cd();
&table[table.size()-1].second);
subc[1]->cd();
subc[2]->cd();
c3c->SaveAs("rhoscanFA.eps");
int iFitBinByBin=table.size();
DoFit(histOutputBBBO,histRhoBBBO,histDataGenO,"bin-by-bin",table);
//=========================== bin-by-bin =================================
//c->cd(1);
//DrawPadProbability(histProbCO);
//c->cd(2);
//DrawPadCorrelations(histRhoBBBO,&table);
c2sq->cd(1);
c2sq->cd(2);
&table[table.size()-1].second);
c2sq->SaveAs("binbybin.eps");
//=========================== iterative ===================================
int iAgoFirstFit=table.size();
for(size_t i=1;i<histRhoAgorep.size();i++) {
int n=nIter[i];
bool isFitted=false;
TString::Format("iterative, N=%d",n),table,n);
isFitted=true;
subc[0]->cd();
TString::Format("iterative N=%d",nIter[i]),0.,
isFitted ? &table[table.size()-1].second : nullptr);
subc[1]->cd();
subc[2]->cd();
c3c->SaveAs(TString::Format("iterative%d.eps",nIter[i]));
}
int iAgoLastFit=table.size();
#ifdef WITH_IDS
int iIDSFirstFit=table.size();
//=========================== IDS ===================================
for(size_t i=2;i<histRhoIDS.size();i++) {
int n=nIter[i];
bool isFitted=false;
TString::Format("IDS, N=%d",n),table,n);
isFitted=true;
subc[0]->cd();
TString::Format("IDS N=%d",nIter[i]),0.,
isFitted ? &table[table.size()-1].second : 0);
subc[1]->cd();
subc[2]->cd();
c3c->SaveAs(TString::Format("ids%d.eps",nIter[i]));
}
int iIDSLastFit=table.size();
#endif
int nfit=table.size();
TH1D *fitChindf=new TH1D("fitChindf",";algorithm;#chi^{2}/NDF",nfit,0,nfit);
TH1D *fitNorm=new TH1D("fitNorm",";algorithm;Landau amplitude [1/GeV]",nfit,0,nfit);
TH1D *fitMu=new TH1D("fitMu",";algorithm;Landau #mu [GeV]",nfit,0,nfit);
TH1D *fitSigma=new TH1D("fitSigma",";algorithm;Landau #sigma [GeV]",nfit,0,nfit);
for(int fit=0;fit<nfit;fit++) {
TF1 *f=table[fit].first;
vector<double> const &r=table[fit].second;
fitChindf->GetXaxis()->SetBinLabel(fit+1,f->GetName());
fitNorm->GetXaxis()->SetBinLabel(fit+1,f->GetName());
fitMu->GetXaxis()->SetBinLabel(fit+1,f->GetName());
fitSigma->GetXaxis()->SetBinLabel(fit+1,f->GetName());
double chi2=r[0];
double ndf=r[1];
fitChindf->SetBinContent(fit+1,chi2/ndf);
fitChindf->SetBinError(fit+1,TMath::Sqrt(2./ndf));
fitNorm->SetBinContent(fit+1,f->GetParameter(0));
fitNorm->SetBinError(fit+1,f->GetParError(0));
fitMu->SetBinContent(fit+1,f->GetParameter(1));
fitMu->SetBinError(fit+1,f->GetParError(1));
fitSigma->SetBinContent(fit+1,f->GetParameter(2));
fitSigma->SetBinError(fit+1,f->GetParError(2));
cout<<"\""<<f->GetName()<<"\","<<r[2]/r[3]<<","<<r[3]
<<","<<TMath::Prob(r[2],r[3])
<<","<<chi2/ndf
<<","<<ndf
<<","<<TMath::Prob(r[0],r[1])
<<","<<r[5];
for(int i=1;i<3;i++) {
cout<<","<<f->GetParameter(i)<<",\"\302\261\","<<f->GetParError(i);
}
cout<<"\n";
}
//=========================== L-curve ==========================
c4->cd(1);
lCurve->SetTitle("L curve;log_{10} L_{x};log_{10} L_{y}");
lCurve->SetLineColor(kRed);
lCurve->Draw("AL");
c4->cd(2);
gPad->Clear();
c4->cd(3);
logTauX->SetTitle(";log_{10} #tau;log_{10} L_{x}");
logTauX->SetLineColor(kBlue);
logTauX->Draw();
c4->cd(4);
logTauY->SetTitle(";log_{10} #tau;log_{10} L_{y}");
logTauY->SetLineColor(kBlue);
logTauY->Draw();
c4->SaveAs("lcurveL.eps");
//========================= rho and L-curve scan ===============
c4->cd(1);
logTauCurvature->SetTitle(";log_{10}(#tau);L curve curvature");
logTauCurvature->SetLineColor(kRed);
logTauCurvature->Draw();
c4->cd(2);
rhoScan->SetTitle(";log_{10}(#tau);average(#rho_{i})");
rhoScan->SetLineColor(kRed);
rhoScan->Draw();
c4->cd(3);
&table[iFitFALCurve].second);
c4->cd(4);
&table[iFitFArho].second);
c4->SaveAs("scanTau.eps");
#ifdef WITH_IDS
#endif
histNiterInversion[0]->GetYaxis()->SetRangeUser(0.3,500.);
histNiterInversion[1]->GetYaxis()->SetRangeUser(-0.1,0.9);
histNiterInversion[2]->GetYaxis()->SetRangeUser(5.6,6.3);
histNiterInversion[3]->GetYaxis()->SetRangeUser(1.6,2.4);
TLine *line=nullptr;
c1->cd();
for(int i=0;i<2;i++) {
gPad->Clear();
gPad->SetLogx();
gPad->SetLogy((i<1));
if(! histNiterInversion[i]) continue;
histNiterInversion[i]->Draw("][");
histNiterFALCurve[i]->Draw("SAME ][");
histNiterFArho[i]->Draw("SAME ][");
histNiterBinByBin[i]->Draw("SAME ][");
graphNiterAgo[i]->Draw("LP");
graphNiterAgoBay[i]->SetMarkerStyle(20);
graphNiterAgoBay[i]->Draw("P");
#ifdef WITH_IDS
graphNiterIDS[i]->Draw("LP");
#endif
if(i==1) {
legend=new TLegend(0.48,0.28,0.87,0.63);
} else {
legend=new TLegend(0.45,0.5,0.88,0.88);
}
legend->SetBorderSize(0);
legend->SetFillStyle(1001);
legend->SetFillColor(kWhite);
legend->SetTextSize(kLegendFontSize*0.7);
legend->AddEntry( histNiterInversion[0],"inversion","l");
legend->AddEntry( histNiterFALCurve[0],"Tikhonov L-curve","l");
legend->AddEntry( histNiterFArho[0],"Tikhonov global cor.","l");
legend->AddEntry( histNiterBinByBin[0],"bin-by-bin","l");
legend->AddEntry( graphNiterAgoBay[0],"\"Bayesian\"","p");
legend->AddEntry( graphNiterAgo[0],"iterative","p");
#ifdef WITH_IDS
legend->AddEntry( graphNiterIDS[0],"IDS","p");
#endif
legend->Draw();
c1->SaveAs(TString::Format("niter%d.eps",i));
}
//c4->cd(1);
//DrawPadCorrelations(histRhoFALCurveO,&table);
c2sq->cd(1);
&table[iFitFALCurve].second,table[iFitFALCurve].first);
c2sq->cd(2);
gPad->SetLogx();
gPad->SetLogy(false);
histNiterInversion[3]->DrawClone("E2");
histNiterInversion[3]->SetFillStyle(0);
histNiterInversion[3]->Draw("SAME HIST ][");
histNiterFALCurve[3]->DrawClone("SAME E2");
histNiterFALCurve[3]->SetFillStyle(0);
histNiterFALCurve[3]->Draw("SAME HIST ][");
histNiterFArho[3]->DrawClone("SAME E2");
histNiterFArho[3]->SetFillStyle(0);
histNiterFArho[3]->Draw("SAME HIST ][");
histNiterBinByBin[3]->DrawClone("SAME E2");
histNiterBinByBin[3]->SetFillStyle(0);
histNiterBinByBin[3]->Draw("SAME HIST ][");
double yTrue=1.8;
line=new TLine(histNiterInversion[3]->GetXaxis()->GetXmin(),
histNiterInversion[3]->GetXaxis()->GetXmax(),
line->Draw();
graphNiterAgo[3]->Draw("LP");
graphNiterAgoBay[3]->SetMarkerStyle(20);
graphNiterAgoBay[3]->Draw("P");
#ifdef WITH_IDS
graphNiterIDS[3]->Draw("LP");
#endif
legend=new TLegend(0.55,0.53,0.95,0.97);
legend->SetBorderSize(0);
legend->SetFillStyle(1001);
legend->SetFillColor(kWhite);
legend->SetTextSize(kLegendFontSize);
legend->AddEntry( line,"truth","l");
legend->AddEntry( histNiterInversion[3],"inversion","l");
legend->AddEntry( histNiterFALCurve[3],"Tikhonov L-curve","l");
legend->AddEntry( histNiterFArho[3],"Tikhonov global cor.","l");
legend->AddEntry( histNiterBinByBin[3],"bin-by-bin","l");
legend->AddEntry( graphNiterAgoBay[3],"\"Bayesian\"","p");
legend->AddEntry( graphNiterAgo[3],"iterative","p");
#ifdef WITH_IDS
legend->AddEntry( graphNiterIDS[3],"IDS","p");
#endif
legend->Draw();
c2sq->SaveAs("fitSigma.eps");
outputFile->Write();
delete outputFile;
}
void GetNiterGraphs(int iFirst,int iLast,vector<pair<TF1*,
vector<double> > > const &table,int color,
TGraph *graph[4],int style) {
for(int ifit=iFirst;ifit<iLast;ifit++) {
vector<double> const &r=table[ifit].second;
niter(ifit-iFirst)=r[4];
chi2(ifit-iFirst)=r[2]/r[3];
gcor(ifit-iFirst)=r[5];
TF1 const *f=table[ifit].first;
mean(ifit-iFirst)=f->GetParameter(1);
emean(ifit-iFirst)=f->GetParError(1);
sigma(ifit-iFirst)=f->GetParameter(2);
esigma(ifit-iFirst)=f->GetParError(2);
}
graph[0]=new TGraph(niter,chi2);
graph[1]=new TGraph(niter,gcor);
graph[2]=new TGraphErrors(niter,mean,eniter,emean);
for(int g=0;g<4;g++) {
if(graph[g]) {
graph[g]->SetLineColor(color);
graph[g]->SetMarkerColor(color);
}
}
}
void GetNiterHist(int ifit,vector<pair<TF1*,vector<double> > > const &table,
TH1 *hist[4],int color,int style,int fillStyle) {
vector<double> const &r=table[ifit].second;
TF1 const *f=table[ifit].first;
hist[0]=new TH1D(table[ifit].first->GetName()+TString("_chi2"),
";iteration;unfold-truth #chi^{2}/N_{D.F.}",1,0.2,1500.);
hist[0]->SetBinContent(1,r[2]/r[3]);
hist[1]=new TH1D(table[ifit].first->GetName()+TString("_gcor"),
";iteration;avg(#rho_{i})",1,0.2,1500.);
hist[1]->SetBinContent(1,r[5]);
hist[2]=new TH1D(table[ifit].first->GetName()+TString("_mu"),
";iteration;parameter #mu",1,0.2,1500.);
hist[2]->SetBinContent(1,f->GetParameter(1));
hist[2]->SetBinError(1,f->GetParError(1));
hist[3]=new TH1D(table[ifit].first->GetName()+TString("_sigma"),
";iteration;parameter #sigma",1,0.2,1500.);
hist[3]->SetBinContent(1,f->GetParameter(2));
hist[3]->SetBinError(1,f->GetParError(2));
for(int h=0;h<4;h++) {
if(hist[h]) {
hist[h]->SetLineColor(color);
hist[h]->SetLineStyle(style);
if( hist[h]->GetBinError(1)>0.0) {
hist[h]->SetFillColor(color-10);
}
hist[h]->SetMarkerStyle(0);
}
}
}
void CreateHistogramCopies(TH1 *h[3],TUnfoldBinning const *binning) {
TString baseName(h[0]->GetName());
h[1]=binning->CreateHistogram(baseName+"_axis",kTRUE,&binMap);
h[2]=(TH1 *)h[1]->Clone(baseName+"_binw");
Int_t nMax=binning->GetEndBin()+1;
for(Int_t iSrc=0;iSrc<nMax;iSrc++) {
double c=h[0]->GetBinContent(iSrc)+h[1]->GetBinContent(iDest);
double e=TMath::Hypot(h[0]->GetBinError(iSrc),h[1]->GetBinError(iDest));
h[1]->SetBinContent(iDest,c);
h[1]->SetBinError(iDest,e);
h[2]->SetBinContent(iDest,c);
h[2]->SetBinError(iDest,e);
}
for(int iDest=0;iDest<=h[2]->GetNbinsX()+1;iDest++) {
double c=h[2]->GetBinContent(iDest);
double e=h[2]->GetBinError(iDest);
double bw=binning->GetBinSize(iDest);
/* if(bw!=h[2]->GetBinWidth(iDest)) {
cout<<"bin "<<iDest<<" width="<<bw<<" "<<h[2]->GetBinWidth(iDest)
<<"\n";
} */
if(bw>0.0) {
h[2]->SetBinContent(iDest,c/bw);
h[2]->SetBinError(iDest,e/bw);
} else {
}
}
}
h[1]=nullptr;
h[2]=nullptr;
}
TH2 *AddOverflowXY(TH2 *h,double widthX,double widthY) {
// add overflow bin to X-axis
int nx=h->GetNbinsX();
int ny=h->GetNbinsY();
double *xBins=new double[nx+2];
double *yBins=new double[ny+2];
for(int i=1;i<=nx;i++) {
xBins[i-1]=h->GetXaxis()->GetBinLowEdge(i);
}
xBins[nx]=h->GetXaxis()->GetBinUpEdge(nx);
for(int i=1;i<=ny;i++) {
yBins[i-1]=h->GetYaxis()->GetBinLowEdge(i);
}
yBins[ny]=h->GetYaxis()->GetBinUpEdge(ny);
TString name(h->GetName());
name+="U";
TH2 *r=new TH2D(name,h->GetTitle(),nx+1,xBins,ny+1,yBins);
for(int ix=0;ix<=nx+1;ix++) {
for(int iy=0;iy<=ny+1;iy++) {
r->SetBinContent(ix,iy,h->GetBinContent(ix,iy));
r->SetBinError(ix,iy,h->GetBinError(ix,iy));
}
}
delete [] yBins;
delete [] xBins;
return r;
}
// add overflow bin to X-axis
int nx=h->GetNbinsX();
double *xBins=new double[nx+2];
for(int i=1;i<=nx;i++) {
xBins[i-1]=h->GetXaxis()->GetBinLowEdge(i);
}
xBins[nx]=h->GetXaxis()->GetBinUpEdge(nx);
TString name(h->GetName());
name+="U";
TH1 *r=new TH1D(name,h->GetTitle(),nx+1,xBins);
for(int ix=0;ix<=nx+1;ix++) {
r->SetBinContent(ix,h->GetBinContent(ix));
r->SetBinError(ix,h->GetBinError(ix));
}
delete [] xBins;
return r;
}
void DrawOverflowX(TH1 *h,double posy) {
double x1=h->GetXaxis()->GetBinLowEdge(h->GetNbinsX());
double x2=h->GetXaxis()->GetBinUpEdge(h->GetNbinsX());
double y0=h->GetYaxis()->GetBinLowEdge(1);
double y2=h->GetYaxis()->GetBinUpEdge(h->GetNbinsY());;
if(h->GetDimension()==1) {
y0=h->GetMinimum();
y2=h->GetMaximum();
}
double w1=-0.3;
TText *textX=new TText((1.+w1)*x2-w1*x1,(1.-posy)*y0+posy*y2,"Overflow bin");
textX->SetNDC(kFALSE);
textX->SetTextSize(0.05);
textX->SetTextAngle(90.);
textX->Draw();
lineX->Draw();
}
void DrawOverflowY(TH1 *h,double posx) {
double x0=h->GetXaxis()->GetBinLowEdge(1);
double x2=h->GetXaxis()->GetBinUpEdge(h->GetNbinsX());
double y1=h->GetYaxis()->GetBinLowEdge(h->GetNbinsY());;
double y2=h->GetYaxis()->GetBinUpEdge(h->GetNbinsY());;
double w1=-0.3;
TText *textY=new TText((1.-posx)*x0+posx*x2,(1.+w1)*y1-w1*y2,"Overflow bin");
textY->SetNDC(kFALSE);
textY->SetTextSize(0.05);
textY->Draw();
TLine *lineY=new TLine(x0,y1,x2,y1);
lineY->Draw();
}
h->Draw("COLZ");
h->SetTitle("migration probabilities;P_{T}(gen) [GeV];P_{T}(rec) [GeV]");
}
h->SetTitle("efficiency;P_{T}(gen) [GeV];#epsilon");
h->SetLineColor(kBlue);
h->SetMinimum(0.75);
h->SetMaximum(1.0);
h->Draw();
TLegend *legEfficiency=new TLegend(0.3,0.58,0.6,0.75);
legEfficiency->SetBorderSize(0);
legEfficiency->SetFillStyle(0);
legEfficiency->AddEntry(h,"reconstruction","l");
legEfficiency->AddEntry((TObject*)nullptr," efficiency","");
legEfficiency->Draw();
}
//gPad->SetLogy(kTRUE);
double amax=0.0;
for(int i=1;i<=histMcRec->GetNbinsX();i++) {
amax=TMath::Max(amax,histMcRec->GetBinContent(i)
+5.0*histMcRec->GetBinError(i));
amax=TMath::Max(amax,histDataRec->GetBinContent(i)
+2.0*histDataRec->GetBinError(i));
}
histMcRec->SetTitle("Reconstructed;P_{T}(rec);Nevent / GeV");
histMcRec->Draw("HIST");
histMcRec->SetLineColor(kBlue);
histMcRec->SetMinimum(1.0);
histMcRec->SetMaximum(amax);
//histMcbgrRec->SetFillMode(1);
histMcbgrRec->SetLineColor(kBlue-6);
histMcbgrRec->SetFillColor(kBlue-10);
histMcbgrRec->Draw("SAME HIST");
TH1 * histFoldBack=nullptr;
histMcRec->Clone(histDataUnfold->GetName()+TString("_folded"));
int nrec=histFoldBack->GetNbinsX();
if((nrec==histProbability->GetNbinsY())&&
(nrec==histMcbgrRec->GetNbinsX())&&
(nrec==histDataRec->GetNbinsX())
) {
for(int ix=1;ix<=nrec;ix++) {
double sum=0.0;
double sume2=0.0;
for(int iy=0;iy<=histProbability->GetNbinsX()+1;iy++) {
sum += histDataUnfold->GetBinContent(iy)*
histDataUnfold->GetBinWidth(iy)*
histProbability->GetBinContent(iy,ix);
for(int iy2=0;iy2<=histProbability->GetNbinsX()+1;iy2++) {
sume2 += histDataUnfold->GetBinError(iy)*
histDataUnfold->GetBinWidth(iy)*
histProbability->GetBinContent(iy,ix)*
histDataUnfold->GetBinError(iy2)*
histDataUnfold->GetBinWidth(iy2)*
histProbability->GetBinContent(iy2,ix)*
histRhoij->GetBinContent(iy,iy2);
}
}
sum /= histFoldBack->GetBinWidth(ix);
sum += histMcbgrRec->GetBinContent(ix);
histFoldBack->SetBinContent(ix,sum);
histFoldBack->SetBinError(ix,TMath::Sqrt(sume2)
/histFoldBack->GetBinWidth(ix));
}
} else {
cout<<"can not fold back: "<<nrec
<<" "<<histProbability->GetNbinsY()
<<" "<<histMcbgrRec->GetNbinsX()
<<" "<<histDataRec->GetNbinsX()
<<"\n";
exit(0);
}
histFoldBack->SetLineColor(kBlack);
histFoldBack->SetMarkerStyle(0);
histFoldBack->Draw("SAME HIST");
}
histDataRec->SetLineColor(kRed);
histDataRec->SetMarkerColor(kRed);
histDataRec->Draw("SAME");
TLegend *legRec=new TLegend(0.4,0.5,0.68,0.85);
legRec->SetBorderSize(0);
legRec->SetFillStyle(0);
legRec->SetTextSize(kLegendFontSize);
legRec->AddEntry(histMcRec,"MC total","l");
legRec->AddEntry(histMcbgrRec,"background","f");
int ndf=-kNbinC;
double sumD=0.,sumF=0.,chi2=0.;
for(int i=1;i<=histDataRec->GetNbinsX();i++) {
//cout<<histDataRec->GetBinContent(i)<<" "<<histFoldBack->GetBinContent(i)<<" "<<" w="<<histFoldBack->GetBinWidth(i)<<"\n";
sumD+=histDataRec->GetBinContent(i)*histDataRec->GetBinWidth(i);
sumF+=histFoldBack->GetBinContent(i)*histFoldBack->GetBinWidth(i);
double pull=(histFoldBack->GetBinContent(i)-histDataRec->GetBinContent(i))/histDataRec->GetBinError(i);
ndf+=1;
}
legRec->AddEntry(histDataRec,TString::Format("data N_{evt}=%.0f",sumD),"lp");
legRec->AddEntry(histFoldBack,TString::Format("folded N_{evt}=%.0f",sumF),"l");
legRec->AddEntry((TObject*)nullptr,TString::Format("#chi^{2}=%.1f ndf=%d",chi2,ndf),"");
//exit(0);
} else {
legRec->AddEntry(histDataRec,"data","lp");
}
legRec->Draw();
}
char const *text,double tau,vector<double> const *r,
TF1 *f) {
//gPad->SetLogy(kTRUE);
double amin=0.;
double amax=0.;
for(int i=1;i<=histMcsigGen->GetNbinsX();i++) {
amin=TMath::Min(amin,histDataUnfold->GetBinContent(i)
-1.1*histDataUnfold->GetBinError(i));
amax=TMath::Max(amax,histDataUnfold->GetBinContent(i)
+1.1*histDataUnfold->GetBinError(i));
}
amin=TMath::Min(amin,histMcsigGen->GetBinContent(i)
-histMcsigGen->GetBinError(i));
amin=TMath::Min(amin,histDataGen->GetBinContent(i)
-histDataGen->GetBinError(i));
amax=TMath::Max(amax,histMcsigGen->GetBinContent(i)
+10.*histMcsigGen->GetBinError(i));
amax=TMath::Max(amax,histDataGen->GetBinContent(i)
+2.*histDataGen->GetBinError(i));
}
histMcsigGen->SetMinimum(amin);
histMcsigGen->SetMaximum(amax);
histMcsigGen->SetTitle("Truth;P_{T};Nevent / GeV");
histMcsigGen->SetLineColor(kBlue);
histMcsigGen->Draw("HIST");
histDataGen->SetLineColor(kRed);
histDataGen->SetMarkerColor(kRed);
histDataGen->SetMarkerSize(1.0);
histDataGen->Draw("SAME HIST");
histDataUnfold->SetMarkerStyle(21);
histDataUnfold->SetMarkerSize(0.7);
histDataUnfold->Draw("SAME");
}
if(f) {
f->SetLineStyle(kSolid);
f->Draw("SAME");
}
TLegend *legTruth=new TLegend(0.32,0.65,0.6,0.9);
legTruth->SetBorderSize(0);
legTruth->SetFillStyle(0);
legTruth->SetTextSize(kLegendFontSize);
legTruth->AddEntry(histMcsigGen,"MC","l");
if(!histDataUnfold) legTruth->AddEntry((TObject *)nullptr," Landau(5,2)","");
legTruth->AddEntry(histDataGen,"data","l");
if(!histDataUnfold) legTruth->AddEntry((TObject *)nullptr," Landau(6,1.8)","");
if(text) t=text;
else t=histDataUnfold->GetName();
if(tau>0) {
t+=TString::Format(" #tau=%.2g",tau);
}
legTruth->AddEntry(histDataUnfold,t,"lp");
if(r) {
legTruth->AddEntry((TObject *)nullptr,"test wrt data:","");
legTruth->AddEntry((TObject *)nullptr,TString::Format
("#chi^{2}/%d=%.1f prob=%.3f",
(int)(*r)[3],(*r)[2]/(*r)[3],
TMath::Prob((*r)[2],(*r)[3])),"");
}
}
if(f) {
legTruth->AddEntry(f,"fit","l");
}
legTruth->Draw();
TPad *subpad = new TPad("subpad","",0.35,0.29,0.88,0.68);
subpad->SetFillStyle(0);
subpad->Draw();
subpad->cd();
amin=0.;
amax=0.;
int istart=11;
for(int i=istart;i<=histMcsigGen->GetNbinsX();i++) {
amin=TMath::Min(amin,histMcsigGen->GetBinContent(i)
-histMcsigGen->GetBinError(i));
amin=TMath::Min(amin,histDataGen->GetBinContent(i)
-histDataGen->GetBinError(i));
amin=TMath::Min(amin,histDataUnfold->GetBinContent(i)
-histDataUnfold->GetBinError(i));
amax=TMath::Max(amax,histMcsigGen->GetBinContent(i)
+histMcsigGen->GetBinError(i));
amax=TMath::Max(amax,histDataGen->GetBinContent(i)
+histDataGen->GetBinError(i));
amax=TMath::Max(amax,histDataUnfold->GetBinContent(i)
+histDataUnfold->GetBinError(i));
}
copyMcsigGen->GetXaxis()->SetRangeUser
(copyMcsigGen->GetXaxis()->GetBinLowEdge(istart),
copyMcsigGen->GetXaxis()->GetBinUpEdge(copyMcsigGen->GetNbinsX()-1));
copyMcsigGen->SetTitle(";;");
copyMcsigGen->GetYaxis()->SetRangeUser(amin,amax);
copyMcsigGen->Draw("HIST");
copyDataGen->Draw("SAME HIST");
copyDataUnfold->Draw("SAME");
if(f) {
((TF1 *)f->Clone())->Draw("SAME");
}
}
}
vector<pair<TF1*,vector<double> > > const *table) {
h->SetMinimum(-1.);
h->SetMaximum(1.);
h->SetTitle("correlation coefficients;P_{T}(gen) [GeV];P_{T}(gen) [GeV]");
h->Draw("COLZ");
if(table) {
TLegend *legGCor=new TLegend(0.13,0.6,0.5,0.8);
legGCor->SetBorderSize(0);
legGCor->SetFillStyle(0);
legGCor->SetTextSize(kLegendFontSize);
vector<double> const &r=(*table)[table->size()-1].second;
legGCor->AddEntry((TObject *)nullptr,TString::Format("min(#rho_{ij})=%5.2f",r[6]),"");
legGCor->AddEntry((TObject *)nullptr,TString::Format("max(#rho_{ij})=%5.2f",r[7]),"");
legGCor->AddEntry((TObject *)nullptr,TString::Format("avg(#rho_i)=%5.2f",r[5]),"");
legGCor->Draw();
}
}
TH1 *g_fcnHist=nullptr;
if(flag==0) {
cout<<"fcn flag=0: npar="<<npar<<" gin="<<gin<<" par=[";
for(int i=0;i<npar;i++) {
cout<<" "<<u[i];
}
cout<<"]\n";
}
int n=g_fcnMatrix->GetNrows();
double x0=0,y0=0.;
for(int i=0;i<=n;i++) {
double x1;
if(i<1) x1=g_fcnHist->GetXaxis()->GetBinLowEdge(i+1);
else x1=g_fcnHist->GetXaxis()->GetBinUpEdge(i);
double y1=TMath::LandauI((x1-u[1])/u[2]);
if(i>0) {
double iy=u[0]*u[2]*(y1-y0)/(x1-x0);
dy(i-1)=iy-g_fcnHist->GetBinContent(i);
//cout<<"i="<<i<<" iy="<<iy<<" delta="<< dy(i-1)<<"\n";
}
x0=x1;
y0=y1;
//cout<<"i="<<i<<" y1="<<y1<<" x1="<<x1<<"\n";
}
TVectorD Hdy=(*g_fcnMatrix) * dy;
//Hdy.Print();
f=Hdy*dy;
//exit(0);
}
TFitResultPtr DoFit(TH1 *h,TH2 *rho,TH1 *truth,const char *text,
vector<pair<TF1 *,vector<double> > > &table,int niter) {
TString option="IESN";
cout<<h->GetName()<<"\n";
double gcorAvg=0.;
double rhoMin=0.;
double rhoMax=0.;
if(rho) {
int n=h->GetNbinsX()-1; // overflow is included as extra bin, exclude in fit
//g_fcnMatrix=new TMatrixD(n,n);
for(int i=0;i<n;i++) {
for(int j=0;j<n;j++) {
v(i,j)=rho->GetBinContent(i+1,j+1)*
(h->GetBinError(i+1)*h->GetBinError(j+1));
}
}
TVectorD di(ev.GetEigenValues());
for(int i=0;i<n;i++) {
if(di(i)>0.0) {
d(i,i)=1./di(i);
} else {
cout<<"bad eigenvalue i="<<i<<" di="<<di(i)<<"\n";
exit(0);
}
}
TMatrixD O(ev.GetEigenVectors());
int error=0;
for(int i=0;i<n;i++) {
if(TMath::Abs(test(i,i)-1.0)>1.E-7) {
error++;
}
for(int j=0;j<n;j++) {
if(i==j) continue;
if(TMath::Abs(test(i,j)>1.E-7)) error++;
}
}
// calculate global correlation coefficient (all bins)
rhoMin=1.;
rhoMax=-1.;
for(int i=0;i<=n;i++) {
for(int j=0;j<=n;j++) {
double rho_ij=rho->GetBinContent(i+1,j+1);
v1(i,j)=rho_ij*
(h->GetBinError(i+1)*h->GetBinError(j+1));
if(i!=j) {
}
}
}
TMatrixD d1(n+1,n+1);
TVectorD di1(ev1.GetEigenValues());
for(int i=0;i<=n;i++) {
if(di1(i)>0.0) {
d1(i,i)=1./di1(i);
} else {
cout<<"bad eigenvalue i="<<i<<" di1="<<di1(i)<<"\n";
exit(0);
}
}
TMatrixD O1(ev1.GetEigenVectors());
for(int i=0;i<=n;i++) {
double gcor2=1.-1./(vinv1(i,i)*v1(i,i));
if(gcor2>=0.0) {
} else {
cout<<"bad global correlation "<<i<<" "<<gcor2<<"\n";
}
}
gcorAvg /=(n+1);
/* if(error) {
v.Print();
g_fcnMatrix->Print();
exit(0);
} */
//g_fcnMatrix->Invert();
//from: HFitImpl.cxx
// TVirtualFitter::FCNFunc_t userFcn = 0;
// typedef void (* FCNFunc_t )(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
// userFcn = (TVirtualFitter::GetFitter())->GetFCN();
// (TVirtualFitter::GetFitter())->SetUserFunc(f1);
//...
//fitok = fitter->FitFCN( userFcn );
option += "U";
}
double xmax=h->GetXaxis()->GetBinUpEdge(h->GetNbinsX()-1);
TF1 *landau=new TF1(text,"[0]*TMath::Landau(x,[1],[2],0)",
0.,xmax);
landau->SetParameter(0,6000.);
landau->SetParameter(1,5.);
landau->SetParameter(2,2.);
landau->SetParError(0,10.);
landau->SetParError(1,0.5);
landau->SetParError(2,0.1);
TFitResultPtr s=h->Fit(landau,option,nullptr,0.,xmax);
int np=landau->GetNpar();
fcn(np,nullptr,r[0],landau->GetParameters(),0);
r[1]=h->GetNbinsX()-1-landau->GetNpar();
for(int i=0;i<h->GetNbinsX()-1;i++) {
double di=h->GetBinContent(i+1)-truth->GetBinContent(i+1);
for(int j=0;j<h->GetNbinsX()-1;j++) {
double dj=h->GetBinContent(j+1)-truth->GetBinContent(j+1);
r[2]+=di*dj*(*g_fcnMatrix)(i,j);
}
} else {
double pull=di/h->GetBinError(i+1);
r[2]+=pull*pull;
}
r[3]+=1.0;
}
r[4]=niter;
if(!niter) r[4]=0.25;
r[5]=gcorAvg;
r[6]=rhoMin;
r[7]=rhoMax;
if(rho) {
g_fcnHist=nullptr;
delete g_fcnMatrix;
g_fcnMatrix=nullptr;
}
table.push_back(make_pair(landau,r));
return s;
}
#ifdef WITH_IDS
//===================== interface to IDS unfolding code follows here
// contact Bogdan Malescu to find it
#include "ids_code.cc"
int N_=data->GetNrows();
soustr = new TVectorD(N_);
for( Int_t i=0; i<N_; i++ ){ (*soustr)[i] = 0.; }
}
int N_=data->GetNrows();
delete unfres2IDS_;
}
#endif
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define g(i)
Definition RSha256.hxx:105
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
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
@ kOrange
Definition Rtypes.h:67
@ kBlack
Definition Rtypes.h:65
@ kMagenta
Definition Rtypes.h:66
@ kWhite
Definition Rtypes.h:65
@ kCyan
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
@ kAzure
Definition Rtypes.h:67
@ kSolid
Definition TAttLine.h:52
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
constexpr Int_t kInfo
Definition TError.h:45
Int_t gErrorIgnoreLevel
Error handling routines.
Definition TError.cxx:31
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char mode
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t width
Option_t Option_t style
Option_t Option_t TPoint TPoint const char text
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
float xmax
TMatrixT< Double_t > TMatrixD
Definition TMatrixDfwd.h:23
R__EXTERN TStyle * gStyle
Definition TStyle.h:436
TVectorT< Double_t > TVectorD
Definition TVectorDfwd.h:23
#define gPad
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:38
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:40
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:44
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:45
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:42
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:39
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:41
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:233
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
A TGraphErrors is a TGraph with error bars.
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:693
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
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:9208
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:6701
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:9224
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:351
Service class for 2-D histogram classes.
Definition TH2.h:30
Double_t GetBinContent(Int_t binx, Int_t biny) const override
Definition TH2.h:97
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
Use the TLine constructor to create a simple line.
Definition TLine.h:22
TMatrixDSymEigen.
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:292
The most important graphics class in the ROOT system.
Definition TPad.h:28
Random number generator class based on M.
Definition TRandom3.h:27
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
Basic string class.
Definition TString.h:139
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:2378
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1640
Base class for several text objects.
Definition TText.h:22
Binning schemes for use with the unfolding algorithm TUnfoldDensity.
TH1 * CreateHistogram(const char *histogramName, Bool_t originalAxisBinning=kFALSE, Int_t **binMap=nullptr, const char *histogramTitle=nullptr, const char *axisSteering=nullptr) const
create a THxx histogram capable to hold the bins of this binning node and its children
Double_t GetBinSize(Int_t iBin) const
get N-dimensional bin size
Int_t GetEndBin(void) const
last+1 bin of this node (includes children)
An algorithm to unfold distributions from detector to truth level.
@ kEScanTauRhoAvg
average global correlation coefficient (from TUnfold::GetRhoI())
@ kDensityModeNone
no scale factors, matrix L is similar to unity matrix
@ kEConstraintArea
enforce preservation of the area
Definition TUnfold.h:119
@ kEConstraintNone
use no extra constraint
Definition TUnfold.h:116
ERegMode
choice of regularisation scheme
Definition TUnfold.h:123
@ kRegModeSize
regularise the amplitude of the output distribution
Definition TUnfold.h:129
@ kHistMapOutputHoriz
truth level on x-axis of the response matrix
Definition TUnfold.h:146
static TVirtualFitter * Fitter(TObject *obj, Int_t maxpar=25)
Static function returning a pointer to the current fitter.
TLine * line
const Double_t sigma
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
fit(model, train_loader, val_loader, num_epochs, batch_size, optimizer, criterion, save_best, scheduler)
double landau(double x, double mu, double sigma)
Definition MathFuncs.h:355
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Double_t Prob(Double_t chi2, Int_t ndf)
Computation of the probability for a certain Chi-squared (chi2) and number of degrees of freedom (ndf...
Definition TMath.cxx:637
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
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Double_t Hypot(Double_t x, Double_t y)
Returns sqrt(x*x + y*y)
Definition TMath.cxx:59
Double_t LandauI(Double_t x)
Returns the cumulative (lower tail integral) of the Landau distribution function at point x.
Definition TMath.cxx:2844
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
th1 Draw()
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345

Version 17.6, in parallel to changes in TUnfold

This file is part of TUnfold.

TUnfold is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

TUnfold is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with TUnfold. If not, see http://www.gnu.org/licenses/.

Author
Stefan Schmitt DESY, 14.10.2008

Definition in file testUnfold7c.C.