Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
TUnfoldDensity.cxx
Go to the documentation of this file.
1// Author: Stefan Schmitt, Amnon Harel
2// DESY and CERN, 11/08/11
3
4// Version 17.9, parallel to changes in TUnfold
5//
6// History:
7// Version 17.8, new method GetDXDY()
8// Version 17.7, with bug-fix for curvature regularisation
9// Version 17.6, with updated doxygen comments and bug-fixes in TUnfoldBinning
10// Version 17.5, bug fix in TUnfold also corrects GetEmatrixSysUncorr()
11// Version 17.4, in parallel to changes in TUnfoldBinning
12// Version 17.3, in parallel to changes in TUnfoldBinning
13// Version 17.2, with new options 'N' and 'c' for axis regularisation steering
14// Version 17.1, add scan type RhoSquare, small bug fixes with useAxisBinning
15// Version 17.0, support for density regularisation, complex binning schemes, tau scan
16
17/** \class TUnfoldDensity
18An algorithm to unfold distributions from detector to truth level
19
20TUnfoldDensity is used to decompose a measurement y into several sources x,
21given the measurement uncertainties, background b and a matrix of migrations A.
22The method can be applied to a large number of problems,
23where the measured distribution y is a linear superposition
24of several Monte Carlo shapes. Beyond such a simple template fit,
25TUnfoldDensity has an adjustable regularisation term and also supports an
26optional constraint on the total number of events.
27Background sources can be specified, with a normalisation constant and
28normalisation uncertainty. In addition, variants of the response
29matrix may be specified, these are taken to determine systematic
30uncertainties. Complex, multidimensional arrangements of signal and
31background bins are managed with the help of the class TUnfoldBinning.
32
33If you use this software, please consider the following citation
34<br/>
35<b>S.Schmitt, JINST 7 (2012) T10003 [arXiv:1205.6201]</b>
36<br/>
37Detailed documentation and updates are available on
38http://www.desy.de/~sschmitt
39
40<h3>Brief recipy to use TUnfoldSys:</h3>
41<ul>
42<li>Set up binning schemes for the truth and measured
43distributions. The binning schemes may be coded in the XML language,
44for reading use TUnfoldBinningXML.</li>
45<li>A matrix (truth,reconstructed) is given as a two-dimensional histogram
46 as argument to the constructor of TUnfold</li>
47<li>A vector of measurements is given as one-dimensional histogram using
48 the SetInput() method</li>
49<li>Repeated calls to SubtractBackground() to specify background
50sources</li>
51<li>Repeated calls to AddSysError() to specify systematic uncertainties
52<li>The unfolding is performed
53<ul>
54<li>either once with a fixed parameter tau, method DoUnfold(tau)</li>
55<li>or multiple times in a scan to determine the best chouce of tau,
56 method ScanLCurve()</li>
57<li>or multiple times in a scan to determine the best chouce of tau,
58 method ScanTau()</li>
59</ul>
60<li>Unfolding results are retrieved using various GetXXX() methods
61</ul>
62A detailed documentation of the various GetXXX() methods to control
63systematic uncertainties is given with the method TUnfoldSys.
64
65<h3>Why to use complex binning schemes</h3>
66
67in literature on unfolding, the "standard" test case is a
68one-dimensional distribution without underflow or overflow bins.
69The migration matrix is almost diagonal.
70<br/>
71<b>This "standard" case is rarely realized for real problems.</b>
72<br/>
73Often one has to deal with multi-dimensional distributions.
74In addition, there are underflow and overflow bins
75or other background bins, possibly determined with the help of auxillary
76measurements.
77<br/>
78In TUnfoldDensity, such complex binning schemes are handled with the help
79of the class TUnfoldBinning. For both the measurement and the truth
80there is a tree structure. The tree nodes may correspond to single
81bins (e.g. nuisance parameters) or may hold multi-dimensional distributions.
82<br/>
83For example, the "measurement" tree could have two leaves, one for
84the primary distribution and one for auxillary measurements.
85Similarly, the "truth" tree could have two leaves, one for the
86signal and one for the background.
87Each of the leaves may then have a multi-dimensional distribution.
88<br/>
89The class TUnfoldBinning takes care to map all bins of the
90"measurement" to a one-dimensional vector y.
91Similarly, the "truth" bins are mapped to the vector x.
92
93<h3>How to choose the regularisation settings</h3>
94
95In TUnfoldDensity, two methods are implemented to determine tau**2
96<ol>
97<li>ScanLcurve() locate the tau where the L-curve plot has a "kink"
98this function is implemented in the TUnfold class</li>
99<li>ScanTau() finds the solution such that some variable
100(e.g. global correlation coefficient) is minimized.
101This function is implemented in the TUnfoldDensity class</li>
102</ol>
103Each of the algorithms has its own advantages and disadvantages.
104The algorithm (1) does not work if the input data are too similar to the
105MC prediction. Typical no-go cases of the L-curve scan are:
106<ul>
107<li>the number of measurements is too small (e.g. ny=nx)
108<li>the input data have no statistical fluctuations
109 [identical MC events are used to fill the matrix of migrations
110 and the vector y for a "closure test"]
111</ul>
112The algorithm (2) only works if the variable does have a real minimum
113as a function of tau. If global correlations are minimized, the situation
114is as follows:
115The matrix of migration typically introduces negative correlations.
116The area constraint introduces some positive correlation.
117Regularisation on the "size" introduces no correlation.
118Regularisation on 1st or 2nd derivatives adds positive correlations.
119<br/>
120For these reasons, "size" regularisation does not work well with
121the tau-scan: the higher tau, the smaller rho, but there is no minimum.
122As a result, large values of tau (too strong regularisation) are found.
123In contrast, the tau-scan is expected to work better with 1st or 2nd
124derivative regularisation, because at some point the negative
125correlations from migrations are approximately cancelled by the
126positive correlations from the regularisation conditions.
127<br/>
128whichever algorithm is used, the output has to be checked:
129<ol>
130<li>The L-curve should have approximate L-shape
131and the final choice of tau should not be at the very edge of the
132scanned region</li>
133<li>The scan result should have a well-defined minimum and the
134final choice of tau should sit right in the minimum</li>
135</ol>
136*/
137
138/*
139 This file is part of TUnfold.
140
141 TUnfold is free software: you can redistribute it and/or modify
142 it under the terms of the GNU General Public License as published by
143 the Free Software Foundation, either version 3 of the License, or
144 (at your option) any later version.
145
146 TUnfold is distributed in the hope that it will be useful,
147 but WITHOUT ANY WARRANTY; without even the implied warranty of
148 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
149 GNU General Public License for more details.
150
151 You should have received a copy of the GNU General Public License
152 along with TUnfold. If not, see <http://www.gnu.org/licenses/>.
153*/
154
155#include "TUnfoldDensity.h"
156#include <TGraph.h>
157#include <TMath.h>
158#include <TVectorD.h>
159#include <TObjString.h>
160#include <iostream>
161#include <map>
162
163//#define DEBUG
164
165#ifdef DEBUG
166using std::cout;
167#endif
168
170
172{
173 // clean up
174 if(fOwnedOutputBins) delete fOwnedOutputBins;
175 if(fOwnedInputBins) delete fOwnedInputBins;
176 if(fRegularisationConditions) delete fRegularisationConditions;
177}
178
179////////////////////////////////////////////////////////////////////////
180/// only for use by root streamer or derived classes
181///
183{
184 fConstOutputBins=nullptr;
185 fConstInputBins=nullptr;
186 fOwnedOutputBins=nullptr;
187 fOwnedInputBins=nullptr;
189}
190
191////////////////////////////////////////////////////////////////////////
192/// set up response matrix A, uncorrelated uncertainties of A,
193/// regularisation scheme and binning schemes
194///
195/// \param[in] hist_A matrix that describes the migrations
196/// \param[in] histmap mapping of the histogram axes to the unfolding output
197/// \param[in] regmode (default=kRegModeSize) global regularisation mode
198/// \param[in] constraint (default=kEConstraintArea) type of constraint
199/// \param[in] densityMode (default=kDensityModeBinWidthAndUser)
200/// regularisation scale factors to construct the matrix L
201/// \param[in] outputBins (default=nullptr) binning scheme for truth (unfolding output)
202/// \param[in] inputBins (default=nullptr) binning scheme for measurement (unfolding
203/// input)
204/// \param[in] regularisationDistribution (default=nullptr) selectin of
205/// regularized distribution
206/// \param[in] regularisationAxisSteering (default=nullptr) detailed
207/// regularisation steeringfor selected distribution
208///
209/// The parameters <b>hist_A, histmap, constraint</b> are
210/// explained with the TUnfoldSys constructor.
211/// <br/>
212/// The parameters <b>outputBins,inputBins</b> set the binning
213/// schemes. If these arguments are zero, simple binning schemes are
214/// constructed which correspond to the axes of the histogram
215/// <b>hist_A</b>.
216/// <br/>
217/// The parameters
218/// <b>regmode, densityMode, regularisationDistribution, regularisationAxisSteering</b>
219/// together control how the initial matrix L of regularisation conditions
220/// is constructed. as explained in RegularizeDistribution().
225 const char *regularisationAxisSteering) :
226 TUnfoldSys(hist_A,histmap,kRegModeNone,constraint)
227{
229 // set up binning schemes
231 fOwnedOutputBins = nullptr;
232 TAxis const *genAxis,*detAxis;
234 genAxis=hist_A->GetXaxis();
235 detAxis=hist_A->GetYaxis();
236 } else {
237 genAxis=hist_A->GetYaxis();
238 detAxis=hist_A->GetXaxis();
239 }
240 if(!fConstOutputBins) {
241 // underflow and overflow are included in the binning scheme
242 // They may be used on generator level
244 new TUnfoldBinning(*genAxis,1,1);
246 }
247 // check whether binning scheme is valid
248 if(fConstOutputBins->GetParentNode()) {
249 Error("TUnfoldDensity",
250 "Invalid output binning scheme (node is not the root node)");
251 }
253 fOwnedInputBins = nullptr;
254 if(!fConstInputBins) {
255 // underflow and overflow are not included in the binning scheme
256 // They are still used to count events which have not been reconstructed
258 new TUnfoldBinning(*detAxis,0,0);
260 }
261 if(fConstInputBins->GetParentNode()) {
262 Error("TUnfoldDensity",
263 "Invalid input binning scheme (node is not the root node)");
264 }
265 // check whether binning scheme matches with the histogram
266 // in terms of total number of bins
267 Int_t nOut=genAxis->GetNbins();
268 Int_t nOutMappedT=TMath::Abs(fConstOutputBins->GetTH1xNumberOfBins(kTRUE));
269 Int_t nOutMappedF=TMath::Abs(fConstOutputBins->GetTH1xNumberOfBins
271 if((nOutMappedT!= nOut)&&(nOutMappedF!=nOut)) {
272 Error("TUnfoldDensity",
273 "Output binning incompatible number of bins: axis %d binning scheme %d (%d)",
275 }
276 // check whether binning scheme matches with the histogram
277 Int_t nInput=detAxis->GetNbins();
278 Int_t nInputMappedT=TMath::Abs(fConstInputBins->GetTH1xNumberOfBins(kTRUE));
279 Int_t nInputMappedF=TMath::Abs(fConstInputBins->GetTH1xNumberOfBins
282 Error("TUnfoldDensity",
283 "Input binning incompatible number of bins:axis %d binning scheme %d (%d) ",
285 }
286
287 // report detailed list of excluded bins
288 for (Int_t ix = 0; ix <= nOut+1; ix++) {
289 if(fHistToX[ix]<0) {
290 Info("TUnfold","*NOT* unfolding bin %s",
291 (char const *)GetOutputBinName(ix));
292 }
293 }
294
295 // set up the regularisation here
296 if(regmode !=kRegModeNone) {
300 }
301}
302
303////////////////////////////////////////////////////////////////////////
304/// Get bin name of an outpt bin
305///
306/// \param[in] iBinX bin number
307///
308/// Return value: name of the bin. The name is constructed from the
309/// entries in the binning scheme and includes information about the
310/// bin borders etc.
311
316
317////////////////////////////////////////////////////////////////////////
318/// density correction factor for a given bin
319///
320/// \param[in] densityMode type of factor to calculate
321/// \param[in] iBin bin number
322///
323/// return a multiplicative factor, for scaling the regularisation
324/// conditions from this bin.
325/// <br/>
326/// For densityMode=kDensityModeNone the factor is set to unity.
327/// For densityMode=kDensityModeBinWidth
328/// the factor is set to 1/binArea
329/// where the binArea is the product of the bin widths in all
330/// dimensions. If the width of a bin is zero or can not be
331/// determined, the factor is set to zero.
332/// For densityMode=kDensityModeUser the factor is determined from the
333/// method TUnfoldBinning::GetBinFactor().
334/// For densityMode=kDensityModeBinWidthAndUser, the results of
335/// kDensityModeBinWidth and kDensityModeUser are multiplied.
338{
339 Double_t factor=1.0;
343 if(binSize>0.0) factor /= binSize;
344 else factor=0.0;
345 }
348 factor *= fConstOutputBins->GetBinFactor(iBin);
349 }
350 return factor;
351}
352
353////////////////////////////////////////////////////////////////////////
354/// set up regularisation conditions
355///
356/// \param[in] regmode basic regularisation mode (see class TUnfold)
357/// \param[in] densityMode how to apply bin-wise factors
358/// \param[in] distribution name of the TUnfoldBinning node for which
359/// the regularisation conditions shall be set (zero matches all nodes)
360/// \param[in] axisSteering regularisation fine-tuning
361///
362/// <b>axisSteering</b> is a string with several tokens, separated by
363/// a semicolon: "axisName[options];axisName[options];...".
364/// <dl>
365/// <dt>axisName</dt>
366/// <dd>the name of an axis. The special name * matches all.
367/// So the argument <b>distribution</b> selects one (or all)
368/// distributions. Witin the selected distribution(s), steering options may be
369/// specified for each axis (or for all axes) to define the
370/// regularisation conditions.</dd>
371/// <dt>options</dt>
372/// <dd>one or several character as follows<br/>
373/// u : exclude underflow bin from derivatives along this axis<br/>
374/// o : exclude overflow bin from derivatives along this axis<br/>
375/// U : exclude underflow bin<br/>
376/// O : exclude overflow bin<br/>
377/// b : use bin width for derivative calculation<br/>
378/// B : same as 'b', in addition normalize to average bin width<br/>
379/// N : completely exclude derivatives along this axis<br/>
380/// p : axis is periodic (e.g. azimuthal angle), so
381/// include derivatives built from combinations involving bins at
382/// both ends of the axis "wrap around"
383/// </dd>
384/// </dl>
385/// example: <b>axisSteering</b>="*[UOB]" uses bin widths to calculate
386/// derivatives but underflow/overflow bins are not regularized
395
396////////////////////////////////////////////////////////////////////////
397/// recursively add regularisation conditions for this node and its children
398///
399/// \param[in] binning current node
400/// \param[in] regmode regularisation mode
401/// \param[in] densityMode type of regularisation scaling
402/// \param[in] distribution target distribution(s) name
403/// \param[in] axisSteering steering within the target distribution(s)
416
417////////////////////////////////////////////////////////////////////////
418/// regularize the distribution fof the given node
419///
420/// \param[in] binning current node
421/// \param[in] regmode regularisation mode
422/// \param[in] densityMode type of regularisation scaling
423/// \param[in] axisSteering detailed steering for the axes of the distribution
424
426(const TUnfoldBinning *binning,ERegMode regmode,
428{
429#ifdef DEBUG
430 cout<<"TUnfoldDensity::RegularizeOneDistribution node="
431 <<binning->GetName()<<" "<<regmode<<" "<<densityMode
432 <<" "<<(axisSteering ? axisSteering : "")<<"\n";
433#endif
435 fRegularisationConditions=new TUnfoldBinning("regularisation");
436
439
440 // decode steering
442 binning->DecodeAxisSteering(axisSteering,"uUoObBpN",isOptionGiven);
443 // U implies u
445 // O implies o
447 // B implies b
449 // option N is removed if any other option is on
450 for(Int_t i=0;i<7;i++) {
452 }
453 // option "c" does not work with options UuOo
454 if(isOptionGiven[6] & (isOptionGiven[0]|isOptionGiven[2]) ) {
455 Error("RegularizeOneDistribution",
456 "axis steering %s is not valid",axisSteering);
457 }
458#ifdef DEBUG
459 cout<<" "<<isOptionGiven[0]
460 <<" "<<isOptionGiven[1]
461 <<" "<<isOptionGiven[2]
462 <<" "<<isOptionGiven[3]
463 <<" "<<isOptionGiven[4]
464 <<" "<<isOptionGiven[5]
465 <<" "<<isOptionGiven[6]
466 <<" "<<isOptionGiven[7]
467 <<"\n";
468#endif
469 Info("RegularizeOneDistribution","regularizing %s regMode=%d"
470 " densityMode=%d axisSteering=%s",
471 binning->GetName(),(Int_t) regmode,(Int_t)densityMode,
473 Int_t startBin=binning->GetStartBin();
475 std::vector<Double_t> factor(endBin-startBin);
476 [[maybe_unused]] Int_t nbin=0;
477 for(Int_t bin=startBin;bin<endBin;bin++) {
478 factor[bin-startBin]=GetDensityFactor(densityMode,bin);
479 if(factor[bin-startBin] !=0.0) nbin++;
480 }
481#ifdef DEBUG
482 cout<<"initial number of bins "<<nbin<<"\n";
483#endif
484 Int_t dimension=binning->GetDistributionDimension();
485
486 // decide whether to skip underflow/overflow bins
487 nbin=0;
488 for(Int_t bin=startBin;bin<endBin;bin++) {
491 if(uStatus & isOptionGiven[1]) factor[bin-startBin]=0.;
492 if(oStatus & isOptionGiven[3]) factor[bin-startBin]=0.;
493 if(factor[bin-startBin] !=0.0) nbin++;
494 }
495#ifdef DEBUG
496 cout<<"after underflow/overflow bin removal "<<nbin<<"\n";
497#endif
498 if(regmode==kRegModeSize) {
499 Int_t nRegBins=0;
500 // regularize all bins of the distribution, possibly excluding
501 // underflow/overflow bins
502 for(Int_t bin=startBin;bin<endBin;bin++) {
503 if(factor[bin-startBin]==0.0) continue;
504 if(AddRegularisationCondition(bin,factor[bin-startBin])) {
505 nRegBins++;
506 }
507 }
508 if(nRegBins) {
509 thisRegularisationBinning->AddBinning("size",nRegBins);
510 }
512 for(Int_t direction=0;direction<dimension;direction++) {
513 // for each direction
514 Int_t nRegBins=0;
517#ifdef DEBUG
518 cout<<"skip direction "<<direction<<"\n";
519#endif
520 continue;
521 }
526 isOptionGiven[2] & directionMask) : 1.0;
527 for(Int_t bin=startBin;bin<endBin;bin++) {
528 // check whether bin is excluded
529 if(factor[bin-startBin]==0.0) continue;
530 // for each bin, find the neighbour bins
533 Int_t error=binning->GetBinNeighbours
536 if(error) {
537 Error("RegularizeOneDistribution",
538 "invalid option %s (isPeriodic) for axis %s"
539 " (has underflow or overflow)",axisSteering,
541 }
542 if((regmode==kRegModeDerivative)&&(iNext>=0)) {
543 Double_t f0 = -factor[bin-startBin];
544 Double_t f1 = factor[iNext-startBin];
546 if(distNext>0.0) {
549 } else {
550 f0=0.;
551 f1=0.;
552 }
553 }
554 if((f0==0.0)||(f1==0.0)) continue;
556 nRegBins++;
557#ifdef DEBUG
558 std::cout<<"Added Reg: bin "<<bin<<" "<<f0
559 <<" next: "<<iNext<<" "<<f1<<"\n";
560#endif
561 }
562 } else if((regmode==kRegModeCurvature)&&(iPrev>=0)&&(iNext>=0)) {
563 Double_t f0 = factor[iPrev-startBin];
564 Double_t f1 = -2.*factor[bin-startBin];
565 Double_t f2 = factor[iNext-startBin];
567 if((distPrev<0.)&&(distNext>0.)) {
571 f0 *= f/distPrev;
572 f1 *= f*(0.5/distPrev+0.5/distNext);
573 f2 *= f/distNext;
574 } else {
575 f0=0.;
576 f1=0.;
577 f2=0.;
578 }
579 }
580 if((f0==0.0)||(f1==0.0)||(f2==0.0)) continue;
582 nRegBins++;
583#ifdef DEBUG
584 std::cout<<"Added Reg: prev "<<iPrev<<" "<<f0
585 <<" bin: "<<bin<<" "<<f1
586 <<" next: "<<iNext<<" "<<f2<<"\n";
587#endif
588 }
589 }
590 }
591 if(nRegBins) {
594 name="derivative_";
595 } else if(regmode==kRegModeCurvature) {
596 name="curvature_";
597 }
600 }
601 }
602 }
603#ifdef DEBUG
604 //fLsquared->Print();
605#endif
606}
607
608///////////////////////////////////////////////////////////////////////
609/// retreive unfolding result as a new histogram
610///
611/// \param[in] histogramName name of the histogram
612/// \param[in] histogramTitle (default=nullptr) title of the histogram
613/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
614/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
615/// distribution
616/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
617/// proper binning and axis labels
618///
619/// return value: pointer to a new histogram. If
620/// <b>useAxisBinning</b> is set and if the selected distribution fits
621/// into a root histogram (1,2,3 dimensions) then return a histogram
622/// with the proper binning on each axis. Otherwise, return a 1D
623/// histogram with equidistant binning. If the histogram title is
624/// zero, a title is assigned automatically.
625///
626/// The <b>axisSteering</b> is defines as follows: "axis[mode];axis[mode];..."
627/// where:
628/// <ul>
629/// <li>axis = name of an axis or *</li>
630/// <li>mode = any combination of the letters CUO0123456789
631/// <ul>
632/// <li>C collapse axis into one bin (add up results). If
633/// any of the numbers 0-9 are given in addition, only these bins are added up.
634/// Here bins are counted from zero, whereas in root, bins are counted
635/// from 1. Obviously, this only works for up to 10 bins.</li>
636/// <li>U discarde underflow bin</li>
637/// <li>O discarde overflow bin</li>
638/// </ul>
639/// </ul>
640/// examples: imagine the binning has two axis, named x and y.
641/// <ul>
642/// <li>"*[UO]" exclude underflow and overflow bins for all axis.
643/// So here a TH2 is returned but all undeflow and overflow bins are empty</li>
644/// <li>"x[UOC123]" integrate over the variable x but only using the
645/// bins 1,2,3 and not the underflow and overflow in x.
646/// Here a TH1 is returned, the axis is labelled "y" and
647/// the underflow and overflow (in y) are filled. However only the x-bins
648/// 1,2,3 are used to determine the content.</li>
649/// <li>"x[C];y[UO]" integrate over the variable x, including
650/// underflow and overflow but exclude underflow and overflow in y.
651/// Here a TH1 is returned, the axis is labelled "y". The underflow
652/// and overflow in y are empty.</li>
653/// </ul>
654
656(const char *histogramName,const char *histogramTitle,
657 const char *distributionName,const char *axisSteering,
659{
660 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
661 Int_t *binMap=nullptr;
662 TH1 *r=binning->CreateHistogram
664 if(r) {
666 }
667 if(binMap) {
668 delete [] binMap;
669 }
670 return r;
671}
672
673////////////////////////////////////////////////////////////////////////
674/// retreive bias vector as a new histogram
675///
676/// \param[in] histogramName name of the histogram
677/// \param[in] histogramTitle (default=nullptr) title of the histogram
678/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
679/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
680/// distribution
681/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
682/// proper binning and axis labels
683///
684/// returns a new histogram. See method GetOutput() for a detailed
685/// description of the arguments
687(const char *histogramName,const char *histogramTitle,
688 const char *distributionName,const char *axisSteering,
690{
691 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
692 Int_t *binMap=nullptr;
693 TH1 *r=binning->CreateHistogram
695 if(r) {
697 }
698 if(binMap) delete [] binMap;
699 return r;
700}
701
702////////////////////////////////////////////////////////////////////////
703/// retreive unfolding result folded back as a new histogram
704///
705/// \param[in] histogramName name of the histogram
706/// \param[in] histogramTitle (default=nullptr) title of the histogram
707/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
708/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
709/// distribution
710/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
711/// proper binning and axis labels
712/// \param[in] addBgr (default=false) if true, include the background
713/// contribution (useful for direct comparison to data)
714///
715/// returns a new histogram. See method GetOutput() for a detailed
716/// description of the arguments
718(const char *histogramName,const char *histogramTitle,
719 const char *distributionName,const char *axisSteering,
721{
722 TUnfoldBinning const *binning=fConstInputBins->FindNode(distributionName);
723 Int_t *binMap=nullptr;
724 TH1 *r=binning->CreateHistogram
726 if(r) {
728 if(addBgr) {
730 }
731 }
732 if(binMap) delete [] binMap;
733 return r;
734}
735
736////////////////////////////////////////////////////////////////////////
737/// retreive a background source in a new histogram
738///
739/// \param[in] histogramName name of the histogram
740/// \param[in] bgrSource the background source to retreive
741/// \param[in] histogramTitle (default=nullptr) title of the histogram
742/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
743/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
744/// distribution
745/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
746/// proper binning and axis labels
747/// \param[in] includeError (default=3) type of background errors to
748/// be included (+1 uncorrelated bgr errors, +2 correlated bgr errors)
749///
750/// returns a new histogram. See method GetOutput() for a detailed
751/// description of the arguments
753(const char *histogramName,const char *bgrSource,const char *histogramTitle,
754 const char *distributionName,const char *axisSteering,Bool_t useAxisBinning,
755 Int_t includeError) const
756{
757 TUnfoldBinning const *binning=fConstInputBins->FindNode(distributionName);
758 Int_t *binMap=nullptr;
759 TH1 *r=binning->CreateHistogram
761 if(r) {
763 }
764 if(binMap) delete [] binMap;
765 return r;
766}
767
768////////////////////////////////////////////////////////////////////////
769/// retreive input distribution in a new histogram
770///
771/// \param[in] histogramName name of the histogram
772/// \param[in] histogramTitle (default=nullptr) title of the histogram
773/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
774/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
775/// distribution
776/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
777/// proper binning and axis labels
778///
779/// returns a new histogram. See method GetOutput() for a detailed
780/// description of the arguments
781///
783(const char *histogramName,const char *histogramTitle,
784 const char *distributionName,const char *axisSteering,
786{
787 TUnfoldBinning const *binning=fConstInputBins->FindNode(distributionName);
788 Int_t *binMap=nullptr;
789 TH1 *r=binning->CreateHistogram
791 if(r) {
793 }
794 if(binMap) delete [] binMap;
795 return r;
796}
797
798////////////////////////////////////////////////////////////////////////
799/// retreive global correlation coefficients including all uncertainty sources
800///
801/// \param[in] histogramName name of the histogram
802/// \param[in] histogramTitle (default=nullptr) title of the histogram
803/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
804/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
805/// distribution
806/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
807/// proper binning and axis labels
808/// \param[out] ematInv (default=nullptr) to return the inverse covariance matrix
809///
810/// returns a new histogram. See method GetOutput() for a detailed
811/// description of the arguments. The inverse of the covariance matrix
812/// is stored in a new histogram returned by <b>ematInv</b> if that
813/// pointer is non-zero.
815(const char *histogramName,const char *histogramTitle,
816 const char *distributionName,const char *axisSteering,
818 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
819 Int_t *binMap=nullptr;
820 TH1 *r=binning->CreateHistogram
822 if(r) {
823 TH2 *invEmat=nullptr;
824 if(ematInv) {
825 if(r->GetDimension()==1) {
827 ematName += "_inverseEMAT";
828 Int_t *binMap2D=nullptr;
832 if(binMap2D) delete [] binMap2D;
833 } else {
834 Error("GetRhoItotal",
835 "can not return inverse of error matrix for this binning");
836 }
837 }
839 if(invEmat) {
841 }
842 }
843 if(binMap) delete [] binMap;
844 return r;
845}
846
847////////////////////////////////////////////////////////////////////////
848/// retreive global correlation coefficients including input
849/// (statistical) and background uncertainties
850///
851/// \param[in] histogramName name of the histogram
852/// \param[in] histogramTitle (default=nullptr) title of the histogram
853/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
854/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
855/// distribution
856/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
857/// proper binning and axis labels
858/// \param[out] ematInv (default=nullptr) to return the inverse covariance matrix
859///
860/// returns a new histogram. See method GetOutput() for a detailed
861/// description of the arguments. The inverse of the covariance matrix
862/// is stored in a new histogram returned by <b>ematInv</b> if that
863/// pointer is non-zero.
865(const char *histogramName,const char *histogramTitle,
866 const char *distributionName,const char *axisSteering,
868 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
869 Int_t *binMap=nullptr;
870 TH1 *r=binning->CreateHistogram
872 if(r) {
873 TH2 *invEmat=nullptr;
874 if(ematInv) {
875 if(r->GetDimension()==1) {
877 ematName += "_inverseEMAT";
878 Int_t *binMap2D=nullptr;
882 if(binMap2D) delete [] binMap2D;
883 } else {
884 Error("GetRhoItotal",
885 "can not return inverse of error matrix for this binning");
886 }
887 }
889 if(invEmat) {
891 }
892 }
893 if(binMap) delete [] binMap;
894 return r;
895}
896
897////////////////////////////////////////////////////////////////////////
898/// retreive a correlated systematic 1-sigma shift
899///
900/// \param[in] source identifier of the systematic uncertainty source
901/// \param[in] histogramName name of the histogram
902/// \param[in] histogramTitle (default=nullptr) title of the histogram
903/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
904/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
905/// distribution
906/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
907/// proper binning and axis labels
908///
909/// returns a new histogram. See method GetOutput() for a detailed
910/// description of the arguments
911
913(const char *source,const char *histogramName,
914 const char *histogramTitle,const char *distributionName,
915 const char *axisSteering,Bool_t useAxisBinning) {
916 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
917 Int_t *binMap=nullptr;
918 TH1 *r=binning->CreateHistogram
920 if(r) {
922 delete r;
923 r=nullptr;
924 }
925 }
926 if(binMap) delete [] binMap;
927 return r;
928}
929
930////////////////////////////////////////////////////////////////////////
931/// retreive systematic 1-sigma shift corresponding to a background
932/// scale uncertainty
933///
934/// \param[in] bgrSource identifier of the background
935/// \param[in] histogramName name of the histogram
936/// \param[in] histogramTitle (default=nullptr) title of the histogram
937/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
938/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
939/// distribution
940/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
941/// proper binning and axis labels
942///
943/// returns a new histogram. See method GetOutput() for a detailed
944/// description of the arguments
946(const char *bgrSource,const char *histogramName,
947 const char *histogramTitle,const char *distributionName,
948 const char *axisSteering,Bool_t useAxisBinning) {
949 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
950 Int_t *binMap=nullptr;
951 TH1 *r=binning->CreateHistogram
953 if(r) {
955 delete r;
956 r=nullptr;
957 }
958 }
959 if(binMap) delete [] binMap;
960 return r;
961}
962
963////////////////////////////////////////////////////////////////////////
964////////////////////////////////////////////////////////////////////////
965/// retreive1-sigma shift corresponding to the previously specified uncertainty on tau
966///
967/// \param[in] histogramName name of the histogram
968/// \param[in] histogramTitle (default=nullptr) title of the histogram
969/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
970/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
971/// distribution
972/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
973/// proper binning and axis labels
974///
975/// returns a new histogram. See method GetOutput() for a detailed
976/// description of the arguments
978(const char *histogramName,const char *histogramTitle,
979 const char *distributionName,const char *axisSteering,Bool_t useAxisBinning)
980{
981 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
982 Int_t *binMap=nullptr;
983 TH1 *r=binning->CreateHistogram
985 if(r) {
987 delete r;
988 r=nullptr;
989 }
990 }
991 if(binMap) delete [] binMap;
992 return r;
993}
994
995////////////////////////////////////////////////////////////////////////
996/// retreive correlation coefficients, including all uncertainties
997///
998/// \param[in] histogramName name of the histogram
999/// \param[in] histogramTitle (default=nullptr) title of the histogram
1000/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
1001/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
1002/// distribution
1003/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
1004/// proper binning and axis labels
1005///
1006/// returns a new histogram. See method GetOutput() for a detailed
1007/// description of the arguments
1009(const char *histogramName,const char *histogramTitle,
1010 const char *distributionName,const char *axisSteering,
1012{
1016 if(r) {
1017 for(Int_t i=0;i<=r->GetNbinsX()+1;i++) {
1018 Double_t e_i=r->GetBinContent(i,i);
1019 if(e_i>0.0) e_i=TMath::Sqrt(e_i);
1020 else e_i=0.0;
1021 for(Int_t j=0;j<=r->GetNbinsY()+1;j++) {
1022 if(i==j) continue;
1023 Double_t e_j=r->GetBinContent(j,j);
1024 if(e_j>0.0) e_j=TMath::Sqrt(e_j);
1025 else e_j=0.0;
1026 Double_t e_ij=r->GetBinContent(i,j);
1027 if((e_i>0.0)&&(e_j>0.0)) {
1028 r->SetBinContent(i,j,e_ij/e_i/e_j);
1029 } else {
1030 r->SetBinContent(i,j,0.0);
1031 }
1032 }
1033 }
1034 for(Int_t i=0;i<=r->GetNbinsX()+1;i++) {
1035 if(r->GetBinContent(i,i)>0.0) {
1036 r->SetBinContent(i,i,1.0);
1037 } else {
1038 r->SetBinContent(i,i,0.0);
1039 }
1040 }
1041 }
1042 return r;
1043}
1044
1045////////////////////////////////////////////////////////////////////////
1046/// retreive covaraince contribution from uncorrelated (statistical)
1047/// uncertainties of the response matrix
1048///
1049/// \param[in] histogramName name of the histogram
1050/// \param[in] histogramTitle (default=nullptr) title of the histogram
1051/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
1052/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
1053/// distribution
1054/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
1055/// proper binning and axis labels
1056///
1057/// returns a new histogram. See method GetOutput() for a detailed
1058/// description of the arguments
1060(const char *histogramName,const char *histogramTitle,
1061 const char *distributionName,const char *axisSteering,
1063{
1064 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
1065 Int_t *binMap=nullptr;
1068 if(r) {
1070 }
1071 if(binMap) delete [] binMap;
1072 return r;
1073}
1074
1075////////////////////////////////////////////////////////////////////////
1076/// retreive covariance contribution from uncorrelated background uncertainties
1077///
1078/// \param[in] bgrSource identifier of the background
1079/// \param[in] histogramName name of the histogram
1080/// \param[in] histogramTitle (default=nullptr) title of the histogram
1081/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
1082/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
1083/// distribution
1084/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
1085/// proper binning and axis labels
1086///
1087/// returns a new histogram. See method GetOutput() for a detailed
1088/// description of the arguments
1090(const char *bgrSource,const char *histogramName,
1091 const char *histogramTitle,const char *distributionName,
1093{
1094 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
1095 Int_t *binMap=nullptr;
1098 if(r) {
1100 }
1101 if(binMap) delete [] binMap;
1102 return r;
1103}
1104
1105////////////////////////////////////////////////////////////////////////
1106/// get covariance contribution from the input uncertainties (data
1107/// statistical uncertainties)
1108///
1109/// \param[in] histogramName name of the histogram
1110/// \param[in] histogramTitle (default=nullptr) title of the histogram
1111/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
1112/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
1113/// distribution
1114/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
1115/// proper binning and axis labels
1116///
1117/// returns a new histogram. See method GetOutput() for a detailed
1118/// description of the arguments
1119
1121(const char *histogramName,const char *histogramTitle,
1122 const char *distributionName,const char *axisSteering,
1124{
1125 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
1126 Int_t *binMap=nullptr;
1129 if(r) {
1131 }
1132 if(binMap) delete [] binMap;
1133 return r;
1134}
1135
1136////////////////////////////////////////////////////////////////////////
1137/// get matrix of probabilities in a new histogram
1138///
1139/// \param[in] histogramName name of the histogram
1140/// \param[in] histogramTitle (default=nullptr) title of the histogram
1141/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
1142/// proper binning and axis labels
1143///
1144/// returns a new histogram. if histogramTitle is null, choose a title
1145/// automatically.
1156
1157////////////////////////////////////////////////////////////////////////
1158/// get matrix DX/DY in a new histogram
1159///
1160/// \param[in] histogramName name of the histogram
1161/// \param[in] histogramTitle (default=nullptr) title of the histogram
1162/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
1163/// proper binning and axis labels
1164///
1165/// returns a new histogram. if histogramTitle is null, choose a title
1166/// automatically.
1177
1178////////////////////////////////////////////////////////////////////////
1179/// get covariance matrix including all contributions
1180///
1181/// \param[in] histogramName name of the histogram
1182/// \param[in] histogramTitle (default=nullptr) title of the histogram
1183/// \param[in] distributionName (default=nullptr) identifier of the distribution to be extracted
1184/// \param[in] axisSteering (default=nullptr) detailed steering within the extracted
1185/// distribution
1186/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
1187/// proper binning and axis labels
1188///
1189/// returns a new histogram. See method GetOutput() for a detailed
1190/// description of the arguments
1191
1193(const char *histogramName,const char *histogramTitle,
1194 const char *distributionName,const char *axisSteering,
1196{
1197 TUnfoldBinning const *binning=fConstOutputBins->FindNode(distributionName);
1198 Int_t *binMap=nullptr;
1201 if(r) {
1203 }
1204 if(binMap) delete [] binMap;
1205 return r;
1206}
1207
1208////////////////////////////////////////////////////////////////////////
1209/// access matrix of regularisation conditions in a new histogram
1210///
1211/// \param[in] histogramName name of the histogram
1212/// \param[in] histogramTitle (default=nullptr) title of the histogram
1213/// \param[in] useAxisBinning (default=true) if set to true, try to extract a histogram with
1214/// proper binning and axis labels
1215///
1216/// returns a new histogram. if histogramTitle is null, choose a title
1217/// automatically.
1218
1220(const char *histogramName,const char *histogramTitle,Bool_t useAxisBinning)
1221{
1225 Warning("GetL",
1226 "remove invalid scheme of regularisation conditions %d %d",
1230 }
1232 fRegularisationConditions=new TUnfoldBinning("regularisation",fL->GetNrows());
1233 Warning("GetL","create flat regularisation conditions scheme");
1234 }
1239 return r;
1240}
1241
1242////////////////////////////////////////////////////////////////////////
1243/// get regularisation conditions multiplied by result vector minus bias
1244/// L(x-biasScale*biasVector)
1245///
1246/// \param[in] histogramName name of the histogram
1247/// \param[in] histogramTitle (default=nullptr) title of the histogram
1248///
1249/// returns a new histogram.
1250/// This is a measure of the level of regulartisation required per
1251/// regularisation condition.
1252/// If there are (negative or positive) spikes,
1253/// these regularisation conditions dominate
1254/// over the other regularisation conditions and may introduce
1255/// the largest biases.
1257(const char *histogramName,const char *histogramTitle)
1258{
1264 Warning("GetLxMinusBias",
1265 "remove invalid scheme of regularisation conditions %d %d",
1269 }
1271 fRegularisationConditions=new TUnfoldBinning("regularisation",fL->GetNrows());
1272 Warning("GetLxMinusBias","create flat regularisation conditions scheme");
1273 }
1276 const Int_t *Ldx_rows=Ldx->GetRowIndexArray();
1277 const Double_t *Ldx_data=Ldx->GetMatrixArray();
1278 for(Int_t row=0;row<Ldx->GetNrows();row++) {
1279 if(Ldx_rows[row]<Ldx_rows[row+1]) {
1280 r->SetBinContent(row+1,Ldx_data[Ldx_rows[row]]);
1281 }
1282 }
1283 delete Ldx;
1284 return r;
1285}
1286
1287////////////////////////////////////////////////////////////////////////
1288/// locate a binning node for the input (measured) quantities
1289///
1290/// \param[in] distributionName (default=nullptr) distribution to look
1291/// for. if zero, return the root node
1292///
1293/// returns: pointer to a TUnfoldBinning object or zero if not found
1295(const char *distributionName) const
1296{
1297 // find binning scheme, input bins
1298 // distributionName : the distribution to locate
1299 return fConstInputBins->FindNode(distributionName);
1300}
1301
1302////////////////////////////////////////////////////////////////////////
1303/// locate a binning node for the unfolded (truth level) quantities
1304///
1305/// \param[in] distributionName (default=nullptr) distribution to look
1306/// for. if zero, return the root node
1307///
1308/// returns: pointer to a TUnfoldBinning object or zero if not found
1309
1311(const char *distributionName) const
1312{
1313 // find binning scheme, output bins
1314 // distributionName : the distribution to locate
1315 return fConstOutputBins->FindNode(distributionName);
1316}
1317
1318////////////////////////////////////////////////////////////////////////
1319/// scan a function wrt tau and determine the minimum
1320///
1321/// \param[in] nPoint number of points to be scanned
1322/// \param[in] tauMin smallest tau value to study
1323/// \param[in] tauMax largest tau value to study
1324/// \param[out] scanResult the scanned function wrt log(tau)
1325/// \param[in] mode 1st parameter for the scan function
1326/// \param[in] distribution 2nd parameter for the scan function
1327/// \param[in] projectionMode 3rd parameter for the scan function
1328/// \param[out] lCurvePlot for monitoring, shows the L-curve
1329/// \param[out] logTauXPlot for monitoring, L-curve(X) as a function of log(tau)
1330/// \param[out] logTauYPlot for monitoring, L-curve(Y) as a function of log(tau)
1331///
1332/// Return value: the coordinate number on the curve <b>scanResult</b>
1333/// which corresponds to the minimum
1334/// <br/>
1335/// The function is scanned by repeating the following steps <b>nPoint</b>
1336/// times
1337/// <ol>
1338/// <li>Choose a value of tau</li>
1339/// <li>Perform the unfolding for this choice of tau, DoUnfold(tau)</li>
1340/// <li>Determinethe scan variable GetScanVariable()</li>
1341/// </ol>
1342/// The method GetScanVariable() defines scans of correlation
1343/// coefficients, where <b>mode</b> is chosen from the enum
1344/// EScanTauMode. In addition one may set <b>distribution</b>
1345/// and/or <b>projectionMode</b> to refine the calculation of
1346/// correlations (e.g. restrict the calcuation to the signal
1347/// distribution and/or exclude underflow and overflow bins).
1348/// See the documentation of GetScanVariable() for details.
1349/// Alternative scan variables may be defined by overriding the
1350/// GetScanVariable() method.
1351/// <br>
1352/// Automatic choice of scan range: if (tauMin,tauMax) do not
1353/// correspond to a valid tau range (e.g. tauMin=tauMax=0.0) then
1354/// the tau range is determined automatically. Use with care!
1357 Int_t mode,const char *distribution,const char *axisSteering,
1359{
1360 typedef std::map<Double_t,Double_t> TauScan_t;
1361 typedef std::map<Double_t,std::pair<Double_t,Double_t> > LCurve_t;
1364
1365 //==========================================================
1366 // algorithm:
1367 // (1) do the unfolding for nPoint-1 points
1368 // and store the results in the map
1369 // curve
1370 // (1a) store minimum and maximum tau to curve
1371 // (1b) insert additional points, until nPoint-1 values
1372 // have been calculated
1373 //
1374 // (2) determine the best choice of tau
1375 // do the unfolding for this point
1376 // and store the result in
1377 // curve
1378 // (3) return the result in scanResult
1379
1380 //==========================================================
1381 // (1) do the unfolding for nPoint-1 points
1382 // and store the results in
1383 // curve
1384 // (1a) store minimum and maximum tau to curve
1385
1386 if((tauMin<=0)||(tauMax<=0.0)||(tauMin>=tauMax)) {
1387 // here no range is given, has to be determined automatically
1388 // the maximum tau is determined from the chi**2 values
1389 // observed from unfolding without regulatisation
1390
1391 // first unfolding, without regularisation
1392 DoUnfold(0.0);
1393
1394 // if the number of degrees of freedom is too small, create an error
1395 if(GetNdf()<=0) {
1396 Error("ScanTau","too few input bins, NDF<=nullptr %d",GetNdf());
1397 }
1401 Info("ScanTau","logtau=-Infinity y=%lf X=%lf Y=%lf",y0,X0,Y0);
1402 {
1403 // unfolding guess maximum tau and store it
1405 0.5*(TMath::Log10(GetChi2A()+3.*TMath::Sqrt(GetNdf()+1.0))
1406 -GetLcurveY());
1409 Fatal("ScanTau","problem (missing regularisation?) X=%f Y=%f",
1411 }
1413 curve[logTau]=y;
1414 lcurve[logTau]=std::make_pair(GetLcurveX(),GetLcurveY());
1415 Info("ScanTau","logtau=%lf y=%lf X=%lf Y=%lf",logTau,y,
1417 }
1418 // minimum tau is chosen such that it is less than
1419 // 1% different from the case of no regularisation
1420 // here, several points are inserted as needed
1421 while(((int)curve.size()<nPoint-1)&&
1422 ((TMath::Abs(GetLcurveX()-X0)>0.00432)||
1423 (TMath::Abs(GetLcurveY()-Y0)>0.00432))) {
1424 Double_t logTau=(*curve.begin()).first-0.5;
1427 curve[logTau]=y;
1428 lcurve[logTau]=std::make_pair(GetLcurveX(),GetLcurveY());
1429 Info("ScanTay","logtau=%lf y=%lf X=%lf Y=%lf",logTau,y,
1431 }
1432 } else {
1435 if(nPoint>1) {
1436 // insert maximum tau
1439 curve[logTauMax]=y;
1440 lcurve[logTauMax]=std::make_pair(GetLcurveX(),GetLcurveY());
1441 Info("ScanTau","logtau=%lf y=%lf X=%lf Y=%lf",logTauMax,y,
1443 }
1444 // insert minimum tau
1447 curve[logTauMin]=y;
1448 lcurve[logTauMin]=std::make_pair(GetLcurveX(),GetLcurveY());
1449 Info("ScanTau","logtau=%lf y=%lf X=%lf Y=%lf",logTauMin,y,
1451 }
1452
1453 //==========================================================
1454 // (1b) insert additional points, until nPoint-1 values
1455 // have been calculated
1456 while((int)curve.size()<nPoint-1) {
1457 // insert additional points
1458 // points are inserted such that the largest interval in log(tau)
1459 // is divided into two smaller intervals
1460 // however, there is a penalty term for subdividing intervals
1461 // which are far away from the minimum
1463 i0=curve.begin();
1464 // locate minimum
1465 Double_t logTauYMin=(*i0).first;
1466 Double_t yMin=(*i0).second;
1467 for(;i0!=curve.end();i0++) {
1468 if((*i0).second<yMin) {
1469 yMin=(*i0).second;
1470 logTauYMin=(*i0).first;
1471 }
1472 }
1473 // insert additional points such that large log(tau) intervals
1474 // near the minimum rho are divided into two
1475 i0=curve.begin();
1476 i1=i0;
1477 Double_t distMax=0.0;
1478 Double_t logTau=0.0;
1479 for(i1++;i1!=curve.end();i1++) {
1480 Double_t dist;
1481 // check size of rho interval
1482 dist=TMath::Abs((*i0).first-(*i1).first)
1483 // penalty term if distance from rhoMax is large
1484 +0.25*TMath::Power(0.5*((*i0).first+(*i1).first)-logTauYMin,2.)/
1485 ((*curve.rbegin()).first-(*curve.begin()).first)/nPoint;
1486 if((dist<=0.0)||(dist>distMax)) {
1487 distMax=dist;
1488 logTau=0.5*((*i0).first+(*i1).first);
1489 }
1490 i0=i1;
1491 }
1494 curve[logTau]=y;
1495 lcurve[logTau]=std::make_pair(GetLcurveX(),GetLcurveY());
1496 Info("ScanTau","logtau=%lf y=%lf X=%lf Y=%lf",logTau,y,
1498 }
1499
1500 //==========================================================
1501 // (2) determine the best choice of tau
1502 // do the unfolding for this point
1503 // and store the result in
1504 // curve
1505
1506 Double_t cTmin=0.0;
1507 {
1508 Double_t *cTi=new Double_t[curve.size()];
1509 Double_t *cCi=new Double_t[curve.size()];
1510 Int_t n=0;
1511 for(TauScan_t::const_iterator i=curve.begin();i!=curve.end();i++) {
1512 cTi[n]=(*i).first;
1513 cCi[n]=(*i).second;
1514 n++;
1515 }
1516 // create rho Spline
1517 TSpline3 *splineC=new TSpline3("L curve curvature",cTi,cCi,n);
1518 // find the maximum of the curvature
1519 // if the parameter iskip is non-zero, then iskip points are
1520 // ignored when looking for the largest curvature
1521 // (there are problems with the curvature determined from the first
1522 // few points of splineX,splineY in the algorithm above)
1523 Int_t iskip=0;
1524 if(n>3) iskip=1;
1525 if(n>6) iskip=2;
1527 cTmin=cTi[iskip];
1528 for(Int_t i=iskip;i<n-1-iskip;i++) {
1529 // find minimum on this spline section
1530 // check boundary conditions for x[i+1]
1531 Double_t xMin=cTi[i+1];
1532 Double_t yMin=cCi[i+1];
1533 if(cCi[i]<yMin) {
1534 yMin=cCi[i];
1535 xMin=cTi[i];
1536 }
1537 // find minimum for x[i]<x<x[i+1]
1538 // get spline coefficiencts and solve equation
1539 // derivative(x)==nullptr
1540 Double_t x,y,b,c,d;
1541 splineC->GetCoeff(i,x,y,b,c,d);
1542 // coefficiencts of quadratic equation
1543 Double_t m_p_half=-c/(3.*d);
1544 Double_t q=b/(3.*d);
1546 if(discr>=0.0) {
1547 // solution found
1549 Double_t xx;
1550 if(m_p_half>0.0) {
1551 xx = m_p_half + discr;
1552 } else {
1553 xx = m_p_half - discr;
1554 }
1555 Double_t dx=cTi[i+1]-x;
1556 // check first solution
1557 if((xx>0.0)&&(xx<dx)) {
1558 y=splineC->Eval(x+xx);
1559 if(y<yMin) {
1560 yMin=y;
1561 xMin=x+xx;
1562 }
1563 }
1564 // second solution
1565 if(xx !=0.0) {
1566 xx= q/xx;
1567 } else {
1568 xx=0.0;
1569 }
1570 // check second solution
1571 if((xx>0.0)&&(xx<dx)) {
1572 y=splineC->Eval(x+xx);
1573 if(y<yMin) {
1574 yMin=y;
1575 xMin=x+xx;
1576 }
1577 }
1578 }
1579 // check whether this local minimum is a global minimum
1580 if(yMin<cCmin) {
1581 cCmin=yMin;
1582 cTmin=xMin;
1583 }
1584 }
1585 delete splineC;
1586 delete[] cTi;
1587 delete[] cCi;
1588 }
1591 {
1593 curve[logTauFin]=y;
1594 lcurve[logTauFin]=std::make_pair(GetLcurveX(),GetLcurveY());
1595 Info("ScanTau","Result logtau=%lf y=%lf X=%lf Y=%lf",logTauFin,y,
1597 }
1598 //==========================================================
1599 // (3) return the result in
1600 // scanResult lCurve logTauX logTauY
1601
1602 Int_t bestChoice=-1;
1603 if(!curve.empty()) {
1604 Double_t *y=new Double_t[curve.size()];
1605 Double_t *logT=new Double_t[curve.size()];
1606 int n=0;
1607 for( TauScan_t::const_iterator i=curve.begin();i!=curve.end();i++) {
1608 if(logTauFin==(*i).first) {
1609 bestChoice=n;
1610 }
1611 y[n]=(*i).second;
1612 logT[n]=(*i).first;
1613 n++;
1614 }
1615 if(scanResult) {
1616 TString name;
1617 name = TString::Format("scan(%d,",mode);
1619 name += ",";
1621 name +=")";
1622 (*scanResult)=new TSpline3(name+"%log(tau)",logT,y,n);
1623 }
1624 delete[] y;
1625 delete[] logT;
1626 }
1627 if(!lcurve.empty()) {
1628 Double_t *logT=new Double_t[lcurve.size()];
1629 Double_t *x=new Double_t[lcurve.size()];
1630 Double_t *y=new Double_t[lcurve.size()];
1631 Int_t n=0;
1632 for(LCurve_t::const_iterator i=lcurve.begin();i!=lcurve.end();i++) {
1633 logT[n]=(*i).first;
1634 x[n]=(*i).second.first;
1635 y[n]=(*i).second.second;
1636 //cout<<logT[n]<<" "<< x[n]<<" "<<y[n]<<"\n";
1637 n++;
1638 }
1639 if(lCurvePlot) {
1640 *lCurvePlot=new TGraph(n,x,y);
1641 (*lCurvePlot)->SetTitle("L curve");
1642 }
1643 if(logTauXPlot)
1644 *logTauXPlot=new TSpline3("log(chi**2)%log(tau)",logT,x,n);
1645 if(logTauYPlot)
1646 *logTauYPlot=new TSpline3("log(reg.cond)%log(tau)",logT,y,n);
1647 delete [] y;
1648 delete [] x;
1649 delete [] logT;
1650 }
1651 return bestChoice;
1652}
1653
1654////////////////////////////////////////////////////////////////////////
1655/// calculate the function for ScanTau()
1656///
1657/// \param[in] mode the variable to be calculated
1658/// \param[in] distribution distribution for which the variable
1659/// is to be calculated
1660/// \param[in] axisSteering detailed steering for selecting bins on
1661/// the axes of the distribution (see method GetRhoItotal())
1662///
1663/// return value: the scan result for the given choice of tau (for
1664/// which the unfolding was performed prior to calling this method)
1665/// <br/>
1666/// In ScanTau() the unfolding is repeated for various choices of tau.
1667/// For each tau, after unfolding, GetScanVariable() is called to
1668/// determine the scan result for this choice of tau.
1669/// <br/>
1670/// the following modes are implemented
1671/// <ul>
1672/// <li>kEScanTauRhoAvg : average (stat+bgr) global correlation</li>
1673/// <li>kEScanTauRhoSquaredAvg : average (stat+bgr) global correlation squared</li>
1674/// <li>kEScanTauRhoMax : maximum (stat+bgr) global correlation</li>
1675/// <li>kEScanTauRhoAvgSys : average (stat+bgr+sys) global correlation</li>
1676/// <li>kEScanTauRhoAvgSquaredSys : average (stat+bgr+sys) global correlation squared</li>
1677/// <li>kEScanTauRhoMaxSys : maximum (stat+bgr+sys) global
1678/// correlation</li>
1679/// </ul>
1680
1682(Int_t mode,const char *distribution,const char *axisSteering)
1683{
1684 Double_t r=0.0;
1685 TString name("GetScanVariable(");
1686 name += TString::Format("%d,",mode);
1688 name += ",";
1690 name += ")";
1691 TH1 *rhoi=nullptr;
1698 }
1699 if(rhoi) {
1700 Double_t sum=0.0;
1701 Double_t rhoMax=0.0;
1702 Int_t n=0;
1703 for(Int_t i=0;i<=rhoi->GetNbinsX()+1;i++) {
1704 Double_t c=rhoi->GetBinContent(i);
1705 if(c>=0.) {
1706 if(c>rhoMax) rhoMax=c;
1707 sum += c;
1708 n ++;
1709 }
1710 }
1712 r=sum/n;
1713 } else if((mode==kEScanTauRhoSquareAvg)||
1715 r=sum/n;
1716 } else {
1717 r=rhoMax;
1718 }
1719 // cout<<r<<" "<<GetRhoAvg()<<" "<<GetRhoMax()<<"\n";
1720 delete rhoi;
1721 } else {
1722 Fatal("GetScanVariable","mode %d not implemented",mode);
1723 }
1724 return r;
1725}
#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
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
#define ClassImp(name)
Definition Rtypes.h:374
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 char Point_t Rectangle_t WindowAttributes_t Float_t r
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 child
Option_t Option_t TPoint TPoint const char mode
char name[80]
Definition TGX11.cxx:110
float * q
const_iterator begin() const
const_iterator end() const
Class to manage histogram axis.
Definition TAxis.h:32
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
Service class for 2-D histogram classes.
Definition TH2.h:30
Int_t GetNrows() const
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1040
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1054
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1082
Class to create third splines to interpolate knots Arbitrary conditions can be introduced for first a...
Definition TSpline.h:182
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
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
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
Int_t GetDistributionDimension(void) const
query dimension of this node's distribution
Int_t GetDistributionNumberOfBins(void) const
number of bins in the distribution possibly including under/overflow
Int_t GetEndBin(void) const
last+1 bin of this node (includes children)
static TH2D * CreateHistogramOfMigrations(TUnfoldBinning const *xAxis, TUnfoldBinning const *yAxis, char const *histogramName, Bool_t originalXAxisBinning=kFALSE, Bool_t originalYAxisBinning=kFALSE, char const *histogramTitle=nullptr)
create a TH2D histogram capable to hold the bins of the two input binning schemes on the x and y axes...
void GetBinUnderflowOverflowStatus(Int_t iBin, Int_t *uStatus, Int_t *oStatus) const
return bit maps indicating underflow and overflow status
virtual Double_t GetDistributionAverageBinSize(Int_t axis, Bool_t includeUnderflow, Bool_t includeOverflow) const
get average bin size on the specified axis
void DecodeAxisSteering(const char *axisSteering, const char *options, Int_t *isOptionGiven) const
decode axis steering
TString GetDistributionAxisLabel(Int_t axis) const
get name of an axis
TUnfoldBinning * AddBinning(TUnfoldBinning *binning)
add a TUnfoldBinning as the last child of this node
Int_t GetBinNeighbours(Int_t globalBin, Int_t axis, Int_t *prev, Double_t *distPrev, Int_t *next, Double_t *distNext, Bool_t isPeriodic=kFALSE) const
get neighbour bins along the specified axis
TH2D * CreateErrorMatrixHistogram(const char *histogramName, Bool_t originalAxisBinning, Int_t **binMap=nullptr, const char *histogramTitle=nullptr, const char *axisSteering=nullptr) const
create a TH2D histogram capable to hold a covariance matrix
Int_t GetStartBin(void) const
first bin of this node
TUnfoldBinning const * GetChildNode(void) const
first daughter node
An algorithm to unfold distributions from detector to truth level.
void RegularizeOneDistribution(const TUnfoldBinning *binning, ERegMode regmode, EDensityMode densityMode, const char *axisSteering)
regularize the distribution fof the given node
@ kEScanTauRhoAvg
average global correlation coefficient (from TUnfold::GetRhoI())
@ kEScanTauRhoMax
maximum global correlation coefficient (from TUnfold::GetRhoI())
@ kEScanTauRhoSquareAvgSys
average global correlation coefficient squared (from TUnfoldSys::GetRhoItotal())
@ kEScanTauRhoMaxSys
maximum global correlation coefficient (from TUnfoldSys::GetRhoItotal())
@ kEScanTauRhoSquareAvg
average global correlation coefficient squared (from TUnfold::GetRhoI())
@ kEScanTauRhoAvgSys
average global correlation coefficient (from TUnfoldSys::GetRhoItotal())
TH1 * GetRhoItotal(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE, TH2 **ematInv=nullptr)
retreive global correlation coefficients including all uncertainty sources
TH2 * GetEmatrixSysUncorr(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE)
retreive covaraince contribution from uncorrelated (statistical) uncertainties of the response matrix
Double_t GetDensityFactor(EDensityMode densityMode, Int_t iBin) const
density correction factor for a given bin
TH2 * GetRhoIJtotal(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE)
retreive correlation coefficients, including all uncertainties
TString GetOutputBinName(Int_t iBinX) const override
Get bin name of an outpt bin.
const TUnfoldBinning * fConstOutputBins
binning scheme for the output (truth level)
const TUnfoldBinning * GetOutputBinning(const char *distributionName=nullptr) const
locate a binning node for the unfolded (truth level) quantities
TH2 * GetEmatrixInput(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE)
get covariance contribution from the input uncertainties (data statistical uncertainties)
TUnfoldBinning * fRegularisationConditions
binning scheme for the regularisation conditions
TH2 * GetL(const char *histogramName, const char *histogramTitle=nullptr, Bool_t useAxisBinning=kTRUE)
access matrix of regularisation conditions in a new histogram
TUnfoldBinning * fOwnedOutputBins
pointer to output binning scheme if owned by this class
TH2 * GetEmatrixTotal(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE)
get covariance matrix including all contributions
TUnfoldBinning * fOwnedInputBins
pointer to input binning scheme if owned by this class
TH1 * GetDeltaSysBackgroundScale(const char *bgrSource, const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE)
retreive systematic 1-sigma shift corresponding to a background scale uncertainty
void RegularizeDistribution(ERegMode regmode, EDensityMode densityMode, const char *distribution, const char *axisSteering)
set up regularisation conditions
TH2 * GetProbabilityMatrix(const char *histogramName, const char *histogramTitle=nullptr, Bool_t useAxisBinning=kTRUE) const
get matrix of probabilities in a new histogram
TH1 * GetOutput(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE) const
retreive unfolding result as a new histogram
TH1 * GetFoldedOutput(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE, Bool_t addBgr=kFALSE) const
retreive unfolding result folded back as a new histogram
const TUnfoldBinning * GetInputBinning(const char *distributionName=nullptr) const
locate a binning node for the input (measured) quantities
TH1 * GetInput(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE) const
retreive input distribution in a new histogram
TH1 * GetRhoIstatbgr(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE, TH2 **ematInv=nullptr)
retreive global correlation coefficients including input (statistical) and background uncertainties
TH1 * GetDeltaSysSource(const char *source, const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE)
retreive a correlated systematic 1-sigma shift
TH1 * GetDeltaSysTau(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE)
retreive1-sigma shift corresponding to the previously specified uncertainty on tau
TH2 * GetEmatrixSysBackgroundUncorr(const char *bgrSource, const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE)
retreive covariance contribution from uncorrelated background uncertainties
const TUnfoldBinning * fConstInputBins
binning scheme for the input (detector level)
TH1 * GetLxMinusBias(const char *histogramName, const char *histogramTitle=nullptr)
get regularisation conditions multiplied by result vector minus bias L(x-biasScale*biasVector)
void RegularizeDistributionRecursive(const TUnfoldBinning *binning, ERegMode regmode, EDensityMode densityMode, const char *distribution, const char *axisSteering)
recursively add regularisation conditions for this node and its children
TUnfoldDensity(void)
only for use by root streamer or derived classes
TH1 * GetBias(const char *histogramName, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE) const
retreive bias vector as a new histogram
virtual Int_t ScanTau(Int_t nPoint, Double_t tauMin, Double_t tauMax, TSpline **scanResult, Int_t mode=kEScanTauRhoAvg, const char *distribution=nullptr, const char *projectionMode=nullptr, TGraph **lCurvePlot=nullptr, TSpline **logTauXPlot=nullptr, TSpline **logTauYPlot=nullptr)
scan a function wrt tau and determine the minimum
virtual Double_t GetScanVariable(Int_t mode, const char *distribution, const char *projectionMode)
calculate the function for ScanTau()
EDensityMode
choice of regularisation scale factors to cinstruct the matrix L
@ kDensityModeUser
scale factors from user function in TUnfoldBinning
@ kDensityModeBinWidthAndUser
scale factors from multidimensional bin width and user function
@ kDensityModeBinWidth
scale factors from multidimensional bin width
TH1 * GetBackground(const char *histogramName, const char *bgrSource=nullptr, const char *histogramTitle=nullptr, const char *distributionName=nullptr, const char *projectionMode=nullptr, Bool_t useAxisBinning=kTRUE, Int_t includeError=3) const
retreive a background source in a new histogram
An algorithm to unfold distributions from detector to truth level, with background subtraction and pr...
Definition TUnfoldSys.h:59
void GetEmatrixTotal(TH2 *ematrix, const Int_t *binMap=nullptr)
Get total error matrix, summing up all contributions.
void GetRhoItotal(TH1 *rhoi, const Int_t *binMap=nullptr, TH2 *invEmat=nullptr)
Get global correlatiocn coefficients, summing up all contributions.
void GetBackground(TH1 *bgr, const char *bgrSource=nullptr, const Int_t *binMap=nullptr, Int_t includeError=3, Bool_t clearHist=kTRUE) const
get background into a histogram
Bool_t GetDeltaSysBackgroundScale(TH1 *delta, const char *source, const Int_t *binMap=nullptr)
correlated one-sigma shifts from background normalisation uncertainty
Bool_t GetDeltaSysTau(TH1 *delta, const Int_t *binMap=nullptr)
correlated one-sigma shifts from shifting tau
Bool_t GetDeltaSysSource(TH1 *hist_delta, const char *source, const Int_t *binMap=nullptr)
correlated one-sigma shifts correspinding to a given systematic uncertainty
void GetEmatrixInput(TH2 *ematrix, const Int_t *binMap=nullptr, Bool_t clearEmat=kTRUE)
covariance matrix contribution from input measurement uncertainties
void GetEmatrixSysBackgroundUncorr(TH2 *ematrix, const char *source, const Int_t *binMap=nullptr, Bool_t clearEmat=kTRUE)
covariance contribution from background uncorrelated uncertainty
void GetEmatrixSysUncorr(TH2 *ematrix, const Int_t *binMap=nullptr, Bool_t clearEmat=kTRUE)
Covariance contribution from uncorrelated uncertainties of the response matrix.
TArrayI fHistToX
mapping of histogram bins to matrix indices
Definition TUnfold.h:170
virtual Double_t GetLcurveY(void) const
get value on y-axis of L-curve determined in recent unfolding
Definition TUnfold.cxx:3261
TMatrixDSparse * MultiplyMSparseM(const TMatrixDSparse *a, const TMatrixD *b) const
multiply sparse matrix and a non-sparse matrix
Definition TUnfold.cxx:761
virtual Double_t DoUnfold(void)
core unfolding algorithm
Definition TUnfold.cxx:247
TMatrixD * fX0
bias vector x0
Definition TUnfold.h:164
void GetBias(TH1 *bias, const Int_t *binMap=nullptr) const
get bias vector including bias scale
Definition TUnfold.cxx:2936
void GetProbabilityMatrix(TH2 *A, EHistMap histmap) const
get matrix of probabilities
Definition TUnfold.cxx:3011
virtual TString GetOutputBinName(Int_t iBinX) const
Get bin name of an outpt bin.
Definition TUnfold.cxx:1668
Double_t fBiasScale
scale factor for the bias
Definition TUnfold.h:162
void GetFoldedOutput(TH1 *folded, const Int_t *binMap=nullptr) const
get unfolding result on detector level
Definition TUnfold.cxx:2963
Bool_t AddRegularisationCondition(Int_t i0, Double_t f0, Int_t i1=-1, Double_t f1=0., Int_t i2=-1, Double_t f2=0.)
add a row of regularisation conditions to the matrix L
Definition TUnfold.cxx:1919
void GetL(TH2 *l) const
get matrix of regularisation conditions
Definition TUnfold.cxx:3192
const TMatrixD * GetX(void) const
vector of the unfolding result
Definition TUnfold.h:242
EConstraint
type of extra constraint
Definition TUnfold.h:113
virtual Double_t GetLcurveX(void) const
get value on x-axis of L-curve determined in recent unfolding
Definition TUnfold.cxx:3252
Double_t GetRhoI(TH1 *rhoi, const Int_t *binMap=nullptr, TH2 *invEmat=nullptr) const
get global correlation coefficiencts, possibly cumulated over several bins
Definition TUnfold.cxx:3505
ERegMode
choice of regularisation scheme
Definition TUnfold.h:123
@ kRegModeNone
no regularisation, or defined later by RegularizeXXX() methods
Definition TUnfold.h:126
@ kRegModeDerivative
regularize the 1st derivative of the output distribution
Definition TUnfold.h:132
@ kRegModeSize
regularise the amplitude of the output distribution
Definition TUnfold.h:129
@ kRegModeCurvature
regularize the 2nd derivative of the output distribution
Definition TUnfold.h:135
void GetInput(TH1 *inputData, const Int_t *binMap=nullptr) const
Input vector of measurements.
Definition TUnfold.cxx:3070
void GetOutput(TH1 *output, const Int_t *binMap=nullptr) const
get output distribution, possibly cumulated over several bins
Definition TUnfold.cxx:3290
EHistMap
arrangement of axes for the response matrix (TH2 histogram)
Definition TUnfold.h:143
@ kHistMapOutputHoriz
truth level on x-axis of the response matrix
Definition TUnfold.h:146
Double_t GetChi2A(void) const
get χ2A contribution determined in recent unfolding
Definition TUnfold.h:329
const TMatrixDSparse * GetDXDY(void) const
matrix of derivatives dx/dy
Definition TUnfold.h:250
TMatrixDSparse * fL
regularisation conditions L
Definition TUnfold.h:156
Int_t GetNdf(void) const
get number of degrees of freedom determined in recent unfolding
Definition TUnfold.h:339
std::ostream & Info()
Definition hadd.cxx:171
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TF1 * f1
Definition legend1.C:11
Int_t Finite(Double_t x)
Check if it is finite with a mask in order to be consistent in presence of fast math.
Definition TMath.h:774
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
Double_t Log10(Double_t x)
Returns the common (base-10) logarithm of x.
Definition TMath.h:766
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345