104#define TUnfold_VERSION "V17.9" 
  105#define TUnfold_CLASS_VERSION 17 
  308   double GetDF(
void) 
const; 
 
 
bool Bool_t
Boolean (0=false, 1=true) (bool)
int Int_t
Signed integer 4 bytes (int)
double Double_t
Double 8 bytes.
#define ClassDefOverride(name, id)
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 dest
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t src
#define TUnfold_CLASS_VERSION
Array of doubles (64 bits per element).
Array of integers (32 bits per element).
A TGraph is an object made of two arrays X and Y with npoints each.
TH1 is the base class of all histogram classes in ROOT.
Service class for 2-D histogram classes.
Mother of all ROOT objects.
Base class for spline implementation containing the Draw/Paint methods.
An algorithm to unfold distributions from detector to truth level.
TArrayI fHistToX
mapping of histogram bins to matrix indices
double GetDF(void) const
return the effecive number of degrees of freedom See e.g.
TMatrixDSparse * fE
matrix E
TMatrixDSparse * fEinv
matrix E^(-1)
virtual Double_t GetLcurveY(void) const
get value on y-axis of L-curve determined in recent unfolding
TMatrixDSparse * fAx
result x folded back A*x
TMatrixDSparse * MultiplyMSparseM(const TMatrixDSparse *a, const TMatrixD *b) const
multiply sparse matrix and a non-sparse matrix
virtual Double_t DoUnfold(void)
core unfolding algorithm
Double_t fChi2A
chi**2 contribution from (y-Ax)Vyy-1(y-Ax)
TMatrixD * fX0
bias vector x0
double GetSURE(void) const
return Stein's unbiased risk estimator See e.g.
void GetBias(TH1 *bias, const Int_t *binMap=nullptr) const
get bias vector including bias scale
TMatrixDSparse * MultiplyMSparseTranspMSparse(const TMatrixDSparse *a, const TMatrixDSparse *b) const
multiply a transposed Sparse matrix with another Sparse matrix
TMatrixDSparse * MultiplyMSparseMSparseTranspVector(const TMatrixDSparse *m1, const TMatrixDSparse *m2, const TMatrixTBase< Double_t > *v) const
calculate a sparse matrix product M1*V*M2T where the diagonal matrix V is given by a vector
TMatrixDSparse * CreateSparseMatrix(Int_t nrow, Int_t ncol, Int_t nele, Int_t *row, Int_t *col, Double_t *data) const
create a sparse matrix, given the nonzero elements
Int_t RegularizeSize(int bin, Double_t scale=1.0)
add a regularisation condition on the magnitude of a truth bin
Double_t fEpsMatrix
machine accuracy used to determine matrix rank after eigenvalue analysis
const TMatrixDSparse * GetDXDAM(int i) const
matrix contributions of the derivative dx/dA
void GetProbabilityMatrix(TH2 *A, EHistMap histmap) const
get matrix of probabilities
Double_t GetChi2L(void) const
get χ2L contribution determined in recent unfolding
TMatrixDSparse * fVxx
covariance matrix Vxx
Int_t GetNy(void) const
returns the number of measurement bins
Double_t GetRhoMax(void) const
get maximum global correlation determined in recent unfolding
virtual TString GetOutputBinName(Int_t iBinX) const
Get bin name of an outpt bin.
Double_t fBiasScale
scale factor for the bias
virtual Int_t ScanSURE(Int_t nPoint, Double_t tauMin, Double_t tauMax, TGraph **logTauSURE=nullptr, TGraph **df_chi2A=nullptr, TGraph **lCurve=nullptr)
minimize Stein's unbiased risk estimator "SURE" using successive calls to DoUnfold at various tau.
virtual Int_t ScanLcurve(Int_t nPoint, Double_t tauMin, Double_t tauMax, TGraph **lCurve, TSpline **logTauX=nullptr, TSpline **logTauY=nullptr, TSpline **logTauCurvature=nullptr)
scan the L curve, determine tau and unfold at the final value of tau
const TMatrixDSparse * GetDXDtauSquared(void) const
vector of derivative dx/dtauSquared, using internal bin counting
Double_t fRhoAvg
average global correlation coefficient
TMatrixDSparse * fDXDtauSquared
derivative of the result wrt tau squared
static void DeleteMatrix(TMatrixD **m)
delete matrix and invalidate pointer
Int_t GetBinFromRow(int ix) const
converts matrix row to truth histogram bin number
void ClearHistogram(TH1 *h, Double_t x=0.) const
Initialize bin contents and bin errors for a given histogram.
Int_t RegularizeDerivative(int left_bin, int right_bin, Double_t scale=1.0)
add a regularisation condition on the difference of two truth bin
Int_t GetNx(void) const
returns internal number of output (truth) matrix rows
const TMatrixDSparse * GetDXDAZ(int i) const
vector contributions of the derivative dx/dA
static const char * GetTUnfoldVersion(void)
return a string describing the TUnfold version
void SetConstraint(EConstraint constraint)
set type of area constraint
void GetFoldedOutput(TH1 *folded, const Int_t *binMap=nullptr) const
get unfolding result on detector level
Int_t RegularizeBins(int start, int step, int nbin, ERegMode regmode)
add regularisation conditions for a group of bins
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
Int_t RegularizeCurvature(int left_bin, int center_bin, int right_bin, Double_t scale_left=1.0, Double_t scale_right=1.0)
add a regularisation condition on the curvature of three truth bin
void SetBias(const TH1 *bias)
set bias vector
void GetL(TH2 *l) const
get matrix of regularisation conditions
ERegMode fRegMode
type of regularisation
Int_t GetNr(void) const
get number of regularisation conditions
TMatrixDSparse * fVxxInv
inverse of covariance matrix Vxx-1
const TMatrixD * GetX(void) const
vector of the unfolding result
TMatrixD * fX
unfolding result x
EConstraint
type of extra constraint
@ kEConstraintArea
enforce preservation of the area
@ kEConstraintNone
use no extra constraint
virtual Double_t GetLcurveX(void) const
get value on x-axis of L-curve determined in recent unfolding
Double_t GetRhoI(TH1 *rhoi, const Int_t *binMap=nullptr, TH2 *invEmat=nullptr) const
get global correlation coefficiencts, possibly cumulated over several bins
TMatrixDSparse * fVyy
covariance matrix Vyy corresponding to y
const TMatrixDSparse * GetVyyInv(void) const
inverse of covariance matrix of the data y
Double_t GetEpsMatrix(void) const
get numerical accuracy for Eigenvalue analysis when inverting matrices with rank problems
Int_t fNdf
number of degrees of freedom
TArrayD fSumOverY
truth vector calculated from the non-normalized response matrix
ERegMode
choice of regularisation scheme
@ kRegModeNone
no regularisation, or defined later by RegularizeXXX() methods
@ kRegModeDerivative
regularize the 1st derivative of the output distribution
@ kRegModeSize
regularise the amplitude of the output distribution
@ kRegModeCurvature
regularize the 2nd derivative of the output distribution
@ kRegModeMixed
mixed regularisation pattern
void GetInput(TH1 *inputData, const Int_t *binMap=nullptr) const
Input vector of measurements.
void SetEpsMatrix(Double_t eps)
set numerical accuracy for Eigenvalue analysis when inverting matrices with rank problems
const TMatrixDSparse * GetE(void) const
matrix E, using internal bin counting
TVectorD GetSqrtEvEmatrix(void) const
Int_t GetRowFromBin(int ix) const
converts truth histogram bin number to matrix row
void GetOutput(TH1 *output, const Int_t *binMap=nullptr) const
get output distribution, possibly cumulated over several bins
void GetRhoIJ(TH2 *rhoij, const Int_t *binMap=nullptr) const
get correlation coefficiencts, possibly cumulated over several bins
void ErrorMatrixToHist(TH2 *ematrix, const TMatrixDSparse *emat, const Int_t *binMap, Bool_t doClear) const
add up an error matrix, also respecting the bin mapping
TArrayI fXToHist
mapping of matrix indices to histogram bins
const TMatrixDSparse * GetAx(void) const
vector of folded-back result
TMatrixDSparse * fDXDY
derivative of the result wrt dx/dy
TMatrixD * fY
input (measured) data y
TMatrixDSparse * InvertMSparseSymmPos(const TMatrixDSparse *A, Int_t *rank) const
get the inverse or pseudo-inverse of a positive, sparse matrix
TMatrixDSparse * fVyyInv
inverse of the input covariance matrix Vyy-1
Double_t fLXsquared
chi**2 contribution from (x-s*x0)TLTL(x-s*x0)
TMatrixDSparse * fDXDAM[2]
matrix contribution to the of derivative dx_k/dA_ij
Double_t fTauSquared
regularisation parameter tau squared
Int_t GetNpar(void) const
get number of truth parameters determined in recent unfolding
const TMatrixDSparse * GetEinv(void) const
matrix E-1, using internal bin counting
virtual void ClearResults(void)
reset all results
Double_t fRhoMax
maximum global correlation coefficient
void GetEmatrix(TH2 *ematrix, const Int_t *binMap=nullptr) const
get output covariance matrix, possibly cumulated over several bins
TMatrixDSparse * MultiplyMSparseMSparse(const TMatrixDSparse *a, const TMatrixDSparse *b) const
multiply two sparse matrices
EConstraint fConstraint
type of constraint to use for the unfolding
TUnfold(void)
only for use by root streamer or derived classes
EHistMap
arrangement of axes for the response matrix (TH2 histogram)
@ kHistMapOutputVert
truth level on y-axis of the response matrix
@ kHistMapOutputHoriz
truth level on x-axis of the response matrix
void AddMSparse(TMatrixDSparse *dest, Double_t f, const TMatrixDSparse *src) const
add a sparse matrix, scaled by a factor, to another scaled matrix
const TMatrixDSparse * GetVxx(void) const
covariance matrix of the result
const TMatrixDSparse * GetVxxInv(void) const
inverse of covariance matrix of the result
void GetNormalisationVector(TH1 *s, const Int_t *binMap=nullptr) const
histogram of truth bins, determined from suming over the response matrix
TMatrixDSparse * fDXDAZ[2]
vector contribution to the of derivative dx_k/dA_ij
Double_t GetRhoIFromMatrix(TH1 *rhoi, const TMatrixDSparse *eOrig, const Int_t *binMap, TH2 *invEmat) const
void InitTUnfold(void)
initialize data menbers, for use in constructors
Double_t GetTau(void) const
return regularisation parameter
Double_t GetChi2A(void) const
get χ2A contribution determined in recent unfolding
Int_t RegularizeBins2D(int start_bin, int step1, int nbin1, int step2, int nbin2, ERegMode regmode)
add regularisation conditions for 2d unfolding
const TMatrixDSparse * GetDXDY(void) const
matrix of derivatives dx/dy
void GetLsquared(TH2 *lsquared) const
get matrix of regularisation conditions squared
void GetInputInverseEmatrix(TH2 *ematrix)
get inverse of the measurement's covariance matrix
TMatrixDSparse * fA
response matrix A
Double_t GetRhoAvg(void) const
get average global correlation determined in recent unfolding
TMatrixDSparse * fL
regularisation conditions L
virtual Int_t SetInput(const TH1 *hist_y, Double_t scaleBias=0.0, Double_t oneOverZeroError=0.0, const TH2 *hist_vyy=nullptr, const TH2 *hist_vyy_inv=nullptr)
Define input data for subsequent calls to DoUnfold(tau)
Int_t fIgnoredBins
number of input bins which are dropped because they have error=nullptr
Int_t GetNdf(void) const
get number of degrees of freedom determined in recent unfolding