102 TString msg =
"All histograms must have equal dimension.\n";
107 msg +=
"Please start again!";
109 Fatal(
"Init", msg,
"%s" );
158 TString msg =
"All histograms must have equal dimension.\n";
164 msg +=
"Please start again!";
166 Fatal(
"Init", msg,
"%s" );
277 for (
Int_t i=0; i<
fNdim; i++) CSVM(i,i) = 1/CSV(i);
280 TMatrixD mCinv = (CVort*CSVM)*CUort;
289 for(
int i=0; i<
fNdim; i++){
296 for(
int i=0; i<
fNdim; i++){
297 for(
int j=0; j<
fNdim; j++){
299 vbtmp(i) += QT(i,j)*vb(j)/BSV(i);
303 mAtmp(i,j) += QT(i,
m)*mA(
m,j)/BSV(i);
339 if (ASV(i)<ASV(0)*eps) sreg = ASV(0)*eps;
341 vdz(i) = sreg/(sreg*sreg + ASV(k)*ASV(k));
342 Z(i,i) = vdz(i)*vdz(i);
348 TMatrixD W = mCinv*Vort*Z*VortT*mCinv;
356 Xtau(i,j) = vxini(i) * vxini(j) * W(i,j);
360 a += mA(
m,i)*mA(
m,j);
362 if(vxini(i) && vxini(j))
363 Xinv(i,j) = a/vxini(i)/vxini(j);
377 Xtau *= 1./scale/scale;
389 Info(
"Unfold",
"Unfolding param: %i",k+1 );
390 Info(
"Unfold",
"Curvature of weight distribution: %f",
GetCurvature( vw, mCurv ) );
394 for(
int i=1; i<=
fNdim; i++){
416 unfcov->
SetTitle(
"Toy covariance matrix");
417 for(
int i=1; i<=
fNdim; i++)
418 for(
int j=1; j<=
fNdim; j++)
432 else L(iPar,iPar) = 0.0;
437 for (
Int_t k=0; k<iPar; k++)
L(iPar,jPar) -=
L(k,iPar)*
L(k,jPar);
438 if (
L(iPar,iPar)!=0.)
L(iPar,jPar) /=
L(iPar,iPar);
439 else L(iPar,jPar) = 0;
452 for (
int i=1; i<=ntoys; i++) {
462 for (
int j=1; j<=
fNdim; j++) {
480 for (
int i=1; i<=ntoys; i++) {
490 for (
int j=1; j<=
fNdim; j++) {
522 unfcov->
SetTitle(
"Toy covariance matrix");
523 for(
int i=1; i<=
fNdim; i++)
524 for(
int j=1; j<=
fNdim; j++)
534 for (
int i=1; i<=ntoys; i++) {
555 for (
int i=1; i<=ntoys; i++) {
678 if (vec2(i) != 0) quot(i) = vec1(i) / vec2(i);
680 if (zero) quot(i) = 0;
681 else quot(i) = vec1(i);
695 if (vec(i) != 0) quotmat(i,j) = mat(i,j) / vec(i);
697 if (zero) quotmat(i,j) = 0;
698 else quotmat(i,j) = mat(i,j);
711 for (
Int_t i=0; i<vec1.
GetNrows(); i++) res(i) = vec1(i) * vec2(i);
720 return vec*(curv*vec);
735 if (
fDdim == 0)
for (
Int_t i=0; i<ndim; i++) tC(i,i) = 1;
736 else if (
fDdim == 1) {
737 for (
Int_t i=0; i<ndim; i++) {
738 if (i < ndim-1) tC(i,i+1) = 1.0;
742 else if (
fDdim == 2) {
743 for (
Int_t i=0; i<ndim; i++) {
744 if (i > 0) tC(i,i-1) = 1.0;
745 if (i < ndim-1) tC(i,i+1) = 1.0;
749 tC(ndim-1,ndim-1) = -1.0;
751 else if (
fDdim == 3) {
752 for (
Int_t i=1; i<ndim-2; i++) {
760 for (
Int_t i=0; i<ndim; i++) {
761 if (i > 0) tC(i,i-1) = -4.0;
762 if (i < ndim-1) tC(i,i+1) = -4.0;
763 if (i > 1) tC(i,i-2) = 1.0;
764 if (i < ndim-2) tC(i,i+2) = 1.0;
768 tC(ndim-1,ndim-1) = 2.0;
770 tC(ndim-2,ndim-1) = -3.0;
772 tC(ndim-1,ndim-2) = -3.0;
774 tC(ndim-2,ndim-2) = 6.0;
776 else if (
fDdim == 5) {
777 for (
Int_t i=2; i < ndim-3; i++) {
786 else if (
fDdim == 6) {
787 for (
Int_t i = 3; i < ndim - 3; i++) {
799 for (
Int_t i=0; i<ndim; i++) tC(i,i) = tC(i,i) + eps;
802 for (
Int_t i=0; i<ndim; i++) {
803 for (
Int_t j=0; j<ndim; j++) {
804 for (
Int_t k=0; k<ndim; k++) {
805 tCurv(i,j) = tCurv(i,j) + tC(k,i)*tC(k,j);
849 if (
TMath::Abs(mat[i][i])/ymax > eps) ipos[nn++] = i;
854 for (
UInt_t in=0; in<
nn; in++)
for (
UInt_t jn=0; jn<
nn; jn++) matwork(in,jn) = 0;
859 matwork[in][in] = mat[ipos[in]][ipos[in]];
861 matwork[in][jn] = mat[ipos[in]][ipos[jn]];
862 matwork[jn][in] = matwork[in][jn];
870 for (
UInt_t i=0; i<
n; i++)
for (
UInt_t j=0; j<
n; j++) mat[i][j] = 0;
874 mat[ipos[in]][ipos[in]] = matwork[in][in];
876 mat[ipos[in]][ipos[jn]] = matwork[in][jn];
877 mat[ipos[jn]][ipos[in]] = mat[ipos[in]][ipos[jn]];
Bool_t fMatToyMode
Internal switch for covariance matrix propagation.
TH2D * fXtau
Distribution of singular values.
SVD Approach to Data Unfolding.
const TVectorD & GetSig()
static void H2Verr(const TH1D *histo, TVectorD &vec)
Fill 1D histogram errors into vector.
Random number generator class based on M.
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Bool_t fToyMode
Toy MC detector response matrix.
static Double_t GetCurvature(const TVectorD &vec, const TMatrixD &curv)
Compute curvature of vector.
TMatrixT< Element > & Transpose(const TMatrixT< Element > &source)
Transpose matrix source.
RooArgList L(const RooAbsArg &v1)
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
TH1D * fSVHist
Distribution of d (for checking regularization)
virtual void SetSeed(ULong_t seed=0)
Set the random generator sequence if seed is 0 (default value) a TUUID is generated and used to fill ...
TH1D * Unfold(Int_t kreg)
Perform the unfolding with regularisation parameter kreg.
Int_t fKReg
Normalize unfolded spectrum to 1.
void FillCurvatureMatrix(TMatrixD &tCurv, TMatrixD &tC) const
TH2D * GetUnfoldCovMatrix(const TH2D *cov, Int_t ntoys, Int_t seed=1)
Determine for given input error matrix covariance matrix of unfolded spectrum from toy simulation giv...
TH2D * GetXinv() const
Returns the computed inverse of the covariance matrix.
static TMatrixD MatDivVec(const TMatrixD &mat, const TVectorD &vec, Int_t zero=0)
Divide matrix entries by vector.
Bool_t fNormalize
Derivative for curvature matrix.
TH2D * GetBCov() const
Returns the covariance matrix.
TSVDUnfold(const TH1D *bdat, const TH1D *bini, const TH1D *xini, const TH2D *Adet)
Alternative constructor User provides data and MC test spectra, as well as detector response matrix...
TH2D * fXinv
Computed regularized covariance matrix.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
TH2D * GetXtau() const
Returns the computed regularized covariance matrix corresponding to total uncertainties on measured s...
Single Value Decomposition class.
TH1D * GetD() const
Returns d vector (for choosing appropriate regularisation)
TH1D * GetSV() const
Returns singular values vector.
Double_t ComputeChiSquared(const TH1D &truspec, const TH1D &unfspec)
Helper routine to compute chi-squared between distributions using the computed inverse of the covaria...
Int_t fDdim
Truth and reconstructed dimensions.
static void RegularisedSymMatInvert(TMatrixDSym &mat, Double_t eps=1e-3)
naive regularised inversion cuts off small elements
virtual void SetBinError(Int_t bin, Double_t error)
See convention for numbering bins in TH1::GetBin.
static void H2V(const TH1D *histo, TVectorD &vec)
Fill 1D histogram into vector.
Element Sum() const
Compute sum of elements.
const TH1D * fBdat
Computed inverse of covariance matrix.
TMatrixT< Double_t > TMatrixD
TH2D * fToymat
Toy MC histogram.
static void V2H(const TVectorD &vec, TH1D &histo)
Fill vector into 1D histogram.
static void H2M(const TH2D *histo, TMatrixD &mat)
Fill 2D histogram into matrix.
TH2D * GetAdetCovMatrix(Int_t ntoys, Int_t seed=1)
Determine covariance matrix of unfolded spectrum from finite statistics in response matrix using pseu...
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...
static TVectorD VecDiv(const TVectorD &vec1, const TVectorD &vec2, Int_t zero=0)
Divide entries of two vectors.
char * Form(const char *fmt,...)
tomato 1-D histogram with a double per channel (see TH1 documentation)}
static TVectorD CompProd(const TVectorD &vec1, const TVectorD &vec2)
Multiply entries of two vectors.
TH1D * fDHist
Regularisation parameter.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
TMatrixTSym< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
Mother of all ROOT objects.
static void M2H(const TMatrixD &mat, TH2D &histo)
Fill 2D histogram into matrix.
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
virtual void SetTitle(const char *title)
Change (i.e.
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
virtual Int_t GetNbinsX() const
virtual Int_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
Double_t Sqrt(Double_t x)
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
virtual ~TSVDUnfold()
Destructor.
virtual Int_t GetNbinsY() const
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
tomato 2-D histogram with a double per channel (see TH1 documentation)}