Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TSVDUnfold.cxx
Go to the documentation of this file.
1// Author: Kerstin Tackmann, Andreas Hoecker, Heiko Lacker
2
3/**********************************************************************************
4 * *
5 * Project: TSVDUnfold - data unfolding based on Singular Value Decomposition *
6 * Package: ROOT *
7 * Class : TSVDUnfold *
8 * *
9 * Description: *
10 * Single class implementation of SVD data unfolding based on: *
11 * A. Hoecker, V. Kartvelishvili, *
12 * "SVD approach to data unfolding" *
13 * NIM A372, 469 (1996) [hep-ph/9509307] *
14 * *
15 * Authors: *
16 * Kerstin Tackmann <Kerstin.Tackmann@cern.ch> - CERN, Switzerland *
17 * Andreas Hoecker <Andreas.Hoecker@cern.ch> - CERN, Switzerland *
18 * Heiko Lacker <lacker@physik.hu-berlin.de> - Humboldt U, Germany *
19 * *
20 * Copyright (c) 2010: *
21 * CERN, Switzerland *
22 * Humboldt University, Germany *
23 * *
24 **********************************************************************************/
25
26////////////////////////////////////////////////////////////////////////////////
27
28/** \class TSVDUnfold
29 \ingroup Hist
30 SVD Approach to Data Unfolding
31 <p>
32 Reference: <a href="http://arXiv.org/abs/hep-ph/9509307">Nucl. Instrum. Meth. A372, 469 (1996) [hep-ph/9509307]</a>
33 <p>
34 TSVDUnfold implements the singular value decomposition based unfolding method (see reference). Currently, the unfolding of one-dimensional histograms is supported, with the same number of bins for the measured and the unfolded spectrum.
35 <p>
36 The unfolding procedure is based on singular value decomposition of the response matrix. The regularisation of the unfolding is implemented via a discrete minimum-curvature condition.
37 <p>
38 Monte Carlo inputs:
39 <ul>
40 <li><tt>xini</tt>: true underlying spectrum (TH1D, n bins)
41 <li><tt>bini</tt>: reconstructed spectrum (TH1D, n bins)
42 <li><tt>Adet</tt>: response matrix (TH2D, nxn bins)
43 </ul>
44 Consider the unfolding of a measured spectrum <tt>bdat</tt> with covariance matrix <tt>Bcov</tt> (if not passed explicitly, a diagonal covariance will be built given the errors of <tt>bdat</tt>). The corresponding spectrum in the Monte Carlo is given by <tt>bini</tt>, with the true underlying spectrum given by <tt>xini</tt>. The detector response is described by <tt>Adet</tt>, with <tt>Adet</tt> filled with events (not probabilities) with the true observable on the y-axis and the reconstructed observable on the x-axis.
45 <p>
46 The measured distribution can be unfolded for any combination of resolution, efficiency and acceptance effects, provided an appropriate definition of <tt>xini</tt> and <tt>Adet</tt>.<br><br>
47 <p>
48 The unfolding can be performed by
49 \code{.cpp}
50 TSVDUnfold *tsvdunf = new TSVDUnfold( bdat, Bcov, bini, xini, Adet );
51 TH1D* unfresult = tsvdunf->Unfold( kreg );
52 \endcode
53 where <tt>kreg</tt> determines the regularisation of the unfolding. In general, overregularisation (too small <tt>kreg</tt>) will bias the unfolded spectrum towards the Monte Carlo input, while underregularisation (too large <tt>kreg</tt>) will lead to large fluctuations in the unfolded spectrum. The optimal regularisation can be determined following guidelines in <a href="http://arXiv.org/abs/hep-ph/9509307">Nucl. Instrum. Meth. A372, 469 (1996) [hep-ph/9509307]</a> using the distribution of the <tt>|d_i|</tt> that can be obtained by <tt>tsvdunf->GetD()</tt> and/or using pseudo-experiments.
54 <p>
55 Covariance matrices on the measured spectrum (for either the total uncertainties or individual sources of uncertainties) can be propagated to covariance matrices using the <tt>GetUnfoldCovMatrix</tt> method, which uses pseudo experiments for the propagation. In addition, <tt>GetAdetCovMatrix</tt> allows for the propagation of the statistical uncertainties on the response matrix using pseudo experiments. The covariance matrix corresponding to <tt>Bcov</tt> is also computed as described in <a href="http://arXiv.org/abs/hep-ph/9509307">Nucl. Instrum. Meth. A372, 469 (1996) [hep-ph/9509307]</a> and can be obtained from <tt>tsvdunf->GetXtau()</tt> and its (regularisation independent) inverse from <tt>tsvdunf->GetXinv()</tt>. The distribution of singular values can be retrieved using <tt>tsvdunf->GetSV()</tt>.
56 <p>
57 See also the tutorial for a toy example.
58*/
59
60#include <iostream>
61
62#include "TSVDUnfold.h"
63#include "TH1D.h"
64#include "TH2D.h"
65#include "TDecompSVD.h"
66#include "TRandom3.h"
67#include "TMath.h"
68
69
70////////////////////////////////////////////////////////////////////////////////
71/// Alternative constructor
72/// User provides data and MC test spectra, as well as detector response matrix, diagonal covariance matrix of measured spectrum built from the uncertainties on measured spectrum
73
74TSVDUnfold::TSVDUnfold( const TH1D *bdat, const TH1D *bini, const TH1D *xini, const TH2D *Adet ) :
75fNdim (0),
76fDdim (2),
77fNormalize (kFALSE),
78fKReg (-1),
79fDHist (nullptr),
80fSVHist (nullptr),
81fXtau (nullptr),
82fXinv (nullptr),
83fBdat (bdat),
84fBini (bini),
85fXini (xini),
86fAdet (Adet),
87fToyhisto (nullptr),
88fToymat (nullptr),
89fToyMode (kFALSE),
90fMatToyMode (kFALSE)
91{
92 if (bdat->GetNbinsX() != bini->GetNbinsX() ||
93 bdat->GetNbinsX() != xini->GetNbinsX() ||
94 bdat->GetNbinsX() != Adet->GetNbinsX() ||
95 bdat->GetNbinsX() != Adet->GetNbinsY()) {
96 TString msg = "All histograms must have equal dimension.\n";
97 msg += Form( " Found: dim(bdat)=%i\n", bdat->GetNbinsX() );
98 msg += Form( " Found: dim(bini)=%i\n", bini->GetNbinsX() );
99 msg += Form( " Found: dim(xini)=%i\n", xini->GetNbinsX() );
100 msg += Form( " Found: dim(Adet)=%i,%i\n", Adet->GetNbinsX(), Adet->GetNbinsY() );
101 msg += "Please start again!";
102
103 Fatal( "Init", msg, "%s" );
104 }
105
106 fBcov = (TH2D*)fAdet->Clone("bcov");
107
108 for(int i=1; i<=fBdat->GetNbinsX(); i++){
109 fBcov->SetBinContent(i, i, fBdat->GetBinError(i)*fBdat->GetBinError(i));
110 for(int j=1; j<=fBdat->GetNbinsX(); j++){
111 if(i==j) continue;
112 fBcov->SetBinContent(i,j,0.);
113 }
114 }
115 // Get the input histos
116 fNdim = bdat->GetNbinsX();
117 fDdim = 2; // This is the derivative used to compute the curvature matrix
118}
119
120
121////////////////////////////////////////////////////////////////////////////////
122/// Default constructor
123/// Initialisation of TSVDUnfold
124/// User provides data and MC test spectra, as well as detector response matrix and the covariance matrix of the measured distribution
125
126TSVDUnfold::TSVDUnfold( const TH1D *bdat, TH2D* Bcov, const TH1D *bini, const TH1D *xini, const TH2D *Adet ) :
127fNdim (0),
128fDdim (2),
129fNormalize (kFALSE),
130fKReg (-1),
131fDHist (nullptr),
132fSVHist (nullptr),
133fXtau (nullptr),
134fXinv (nullptr),
135fBdat (bdat),
136fBcov (Bcov),
137fBini (bini),
138fXini (xini),
139fAdet (Adet),
140fToyhisto (nullptr),
141fToymat (nullptr),
142fToyMode (kFALSE),
143fMatToyMode (kFALSE)
144{
145 if (bdat->GetNbinsX() != bini->GetNbinsX() ||
146 bdat->GetNbinsX() != xini->GetNbinsX() ||
147 bdat->GetNbinsX() != Bcov->GetNbinsX() ||
148 bdat->GetNbinsX() != Bcov->GetNbinsY() ||
149 bdat->GetNbinsX() != Adet->GetNbinsX() ||
150 bdat->GetNbinsX() != Adet->GetNbinsY()) {
151 TString msg = "All histograms must have equal dimension.\n";
152 msg += Form( " Found: dim(bdat)=%i\n", bdat->GetNbinsX() );
153 msg += Form( " Found: dim(Bcov)=%i,%i\n", Bcov->GetNbinsX(), Bcov->GetNbinsY() );
154 msg += Form( " Found: dim(bini)=%i\n", bini->GetNbinsX() );
155 msg += Form( " Found: dim(xini)=%i\n", xini->GetNbinsX() );
156 msg += Form( " Found: dim(Adet)=%i,%i\n", Adet->GetNbinsX(), Adet->GetNbinsY() );
157 msg += "Please start again!";
158
159 Fatal( "Init", msg, "%s" );
160 }
161
162 // Get the input histos
163 fNdim = bdat->GetNbinsX();
164 fDdim = 2; // This is the derivative used to compute the curvature matrix
165}
166
167////////////////////////////////////////////////////////////////////////////////
168/// Copy constructor
169
171: TObject ( other ),
172fNdim (other.fNdim),
173fDdim (other.fDdim),
174fNormalize (other.fNormalize),
175fKReg (other.fKReg),
176fDHist (other.fDHist),
177fSVHist (other.fSVHist),
178fXtau (other.fXtau),
179fXinv (other.fXinv),
180fBdat (other.fBdat),
181fBcov (other.fBcov),
182fBini (other.fBini),
183fXini (other.fXini),
184fAdet (other.fAdet),
185fToyhisto (other.fToyhisto),
186fToymat (other.fToymat),
187fToyMode (other.fToyMode),
188fMatToyMode (other.fMatToyMode)
189{
190}
191
192////////////////////////////////////////////////////////////////////////////////
193/// Destructor
194
196{
197 if(fToyhisto){
198 delete fToyhisto;
199 fToyhisto = nullptr;
200 }
201
202 if(fToymat){
203 delete fToymat;
204 fToymat = nullptr;
205 }
206
207 if(fDHist){
208 delete fDHist;
209 fDHist = nullptr;
210 }
211
212 if(fSVHist){
213 delete fSVHist;
214 fSVHist = nullptr;
215 }
216
217 if(fXtau){
218 delete fXtau;
219 fXtau = nullptr;
220 }
221
222 if(fXinv){
223 delete fXinv;
224 fXinv = nullptr;
225 }
226
227 if(fBcov){
228 delete fBcov;
229 fBcov = nullptr;
230 }
231}
232
233////////////////////////////////////////////////////////////////////////////////
234/// Perform the unfolding with regularisation parameter kreg
235
237{
238 fKReg = kreg;
239
240 // Make the histos
241 if (!fToyMode && !fMatToyMode) InitHistos( );
242
243 // Create vectors and matrices
246
247 Double_t eps = 1e-12;
249
250 // Copy histograms entries into vector
251 if (fToyMode) { H2V( fToyhisto, vb ); H2Verr( fToyhisto, vberr ); }
252 else { H2V( fBdat, vb ); H2Verr( fBdat, vberr ); }
253
254 H2M( fBcov, mB);
255 H2V( fBini, vbini );
256 H2V( fXini, vxini );
257 if (fMatToyMode) H2M( fToymat, mA );
258 else H2M( fAdet, mA );
259
260 // Fill and invert the second derivative matrix
262
263 // Inversion of mC by help of SVD
265 TMatrixD CUort = CSVD.GetU();
266 TMatrixD CVort = CSVD.GetV();
267 TVectorD CSV = CSVD.GetSig();
268
270 for (Int_t i=0; i<fNdim; i++) CSVM(i,i) = 1/CSV(i);
271
272 CUort.Transpose( CUort );
274
275 //Rescale using the data covariance matrix
276 TDecompSVD BSVD( mB );
277 TMatrixD QT = BSVD.GetU();
278 QT.Transpose(QT);
279 TVectorD B2SV = BSVD.GetSig();
281
282 for(int i=0; i<fNdim; i++){
283 BSV(i) = TMath::Sqrt(B2SV(i));
284 }
287 mAtmp *= 0;
288 vbtmp *= 0;
289 for(int i=0; i<fNdim; i++){
290 for(int j=0; j<fNdim; j++){
291 if(BSV(i)){
292 vbtmp(i) += QT(i,j)*vb(j)/BSV(i);
293 }
294 for(int m=0; m<fNdim; m++){
295 if(BSV(i)){
296 mAtmp(i,j) += QT(i,m)*mA(m,j)/BSV(i);
297 }
298 }
299 }
300 }
301 mA = mAtmp;
302 vb = vbtmp;
303
304 // Singular value decomposition and matrix operations
306 TMatrixD Uort = ASVD.GetU();
307 TMatrixD Vort = ASVD.GetV();
308 TVectorD ASV = ASVD.GetSig();
309
310 if (!fToyMode && !fMatToyMode) {
311 V2H(ASV, *fSVHist);
312 }
313
315
316 Uort.Transpose(Uort);
317 TVectorD vd = Uort*vb;
318
319 if (!fToyMode && !fMatToyMode) {
320 V2H(vd, *fDHist);
321 }
322
323 // Damping coefficient
324 Int_t k = GetKReg()-1;
325
326 TVectorD vx(fNdim); // Return variable
327
328 // Damping factors
330 TMatrixD Z(fNdim, fNdim);
331 for (Int_t i=0; i<fNdim; i++) {
332 if (ASV(i)<ASV(0)*eps) sreg = ASV(0)*eps;
333 else sreg = ASV(i);
334 vdz(i) = sreg/(sreg*sreg + ASV(k)*ASV(k));
335 Z(i,i) = vdz(i)*vdz(i);
336 }
337 TVectorD vz = CompProd( vd, vdz );
338
340 VortT.Transpose(VortT);
342
345 Xtau *= 0;
346 Xinv *= 0;
347 for (Int_t i=0; i<fNdim; i++) {
348 for (Int_t j=0; j<fNdim; j++) {
349 Xtau(i,j) = vxini(i) * vxini(j) * W(i,j);
350
351 double a=0;
352 for (Int_t m=0; m<fNdim; m++) {
353 a += mA(m,i)*mA(m,j);
354 }
355 if(vxini(i) && vxini(j))
356 Xinv(i,j) = a/vxini(i)/vxini(j);
357 }
358 }
359
360 // Compute the weights
361 TVectorD vw = Vreg*vz;
362
363 // Rescale by xini
364 vx = CompProd( vw, vxini );
365
366 if(fNormalize){ // Scale result to unit area
367 Double_t scale = vx.Sum();
368 if (scale > 0){
369 vx *= 1.0/scale;
370 Xtau *= 1./scale/scale;
371 Xinv *= scale*scale;
372 }
373 }
374
375 if (!fToyMode && !fMatToyMode) {
376 M2H(Xtau, *fXtau);
377 M2H(Xinv, *fXinv);
378 }
379
380 // Get Curvature and also chi2 in case of MC unfolding
381 if (!fToyMode && !fMatToyMode) {
382 Info( "Unfold", "Unfolding param: %i",k+1 );
383 Info( "Unfold", "Curvature of weight distribution: %f", GetCurvature( vw, mCurv ) );
384 }
385
386 TH1D* h = (TH1D*)fBdat->Clone("unfoldingresult");
387 for(int i=1; i<=fNdim; i++){
388 h->SetBinContent(i,0.);
389 h->SetBinError(i,0.);
390 }
391 V2H( vx, *h );
392
393 return h;
394}
395
396////////////////////////////////////////////////////////////////////////////////
397/// Determine for given input error matrix covariance matrix of unfolded
398/// spectrum from toy simulation given the passed covariance matrix on measured spectrum
399/// "cov" - covariance matrix on the measured spectrum, to be propagated
400/// "ntoys" - number of pseudo experiments used for the propagation
401/// "seed" - seed for pseudo experiments
402/// Note that this covariance matrix will contain effects of forced normalisation if spectrum is normalised to unit area.
403
405{
406 fToyMode = true;
407 TH1D* unfres = nullptr;
408 TH2D* unfcov = (TH2D*)fAdet->Clone("unfcovmat");
409 unfcov->SetTitle("Toy covariance matrix");
410 for(int i=1; i<=fNdim; i++)
411 for(int j=1; j<=fNdim; j++)
412 unfcov->SetBinContent(i,j,0.);
413
414 // Code for generation of toys (taken from RooResult and modified)
415 // Calculate the elements of the upper-triangular matrix L that
416 // gives Lt*L = C, where Lt is the transpose of L (the "square-root method")
417 TMatrixD L(fNdim,fNdim); L *= 0;
418
419 for (Int_t iPar= 0; iPar < fNdim; iPar++) {
420
421 // Calculate the diagonal term first
422 L(iPar,iPar) = cov->GetBinContent(iPar+1,iPar+1);
423 for (Int_t k=0; k<iPar; k++) L(iPar,iPar) -= TMath::Power( L(k,iPar), 2 );
424 if (L(iPar,iPar) > 0.0) L(iPar,iPar) = TMath::Sqrt(L(iPar,iPar));
425 else L(iPar,iPar) = 0.0;
426
427 // ...then the off-diagonal terms
428 for (Int_t jPar=iPar+1; jPar<fNdim; jPar++) {
429 L(iPar,jPar) = cov->GetBinContent(iPar+1,jPar+1);
430 for (Int_t k=0; k<iPar; k++) L(iPar,jPar) -= L(k,iPar)*L(k,jPar);
431 if (L(iPar,iPar)!=0.) L(iPar,jPar) /= L(iPar,iPar);
432 else L(iPar,jPar) = 0;
433 }
434 }
435
436 // Remember it
438 TRandom3 random(seed);
439
440 fToyhisto = (TH1D*)fBdat->Clone("toyhisto");
441 TH1D *toymean = (TH1D*)fBdat->Clone("toymean");
442 for (Int_t j=1; j<=fNdim; j++) toymean->SetBinContent(j,0.);
443
444 // Get the mean of the toys first
445 for (int i=1; i<=ntoys; i++) {
446
447 // create a vector of unit Gaussian variables
449 for (Int_t k= 0; k < fNdim; k++) g(k) = random.Gaus(0.,1.);
450
451 // Multiply this vector by Lt to introduce the appropriate correlations
452 g *= (*Lt);
453
454 // Add the mean value offsets and store the results
455 for (int j=1; j<=fNdim; j++) {
456 fToyhisto->SetBinContent(j,fBdat->GetBinContent(j)+g(j-1));
457 fToyhisto->SetBinError(j,fBdat->GetBinError(j));
458 }
459
460 unfres = Unfold(GetKReg());
461
462 for (Int_t j=1; j<=fNdim; j++) {
463 toymean->SetBinContent(j, toymean->GetBinContent(j) + unfres->GetBinContent(j)/ntoys);
464 }
465 delete unfres;
466 unfres = nullptr;
467 }
468
469 // Reset the random seed
470 random.SetSeed(seed);
471
472 //Now the toys for the covariance matrix
473 for (int i=1; i<=ntoys; i++) {
474
475 // Create a vector of unit Gaussian variables
477 for (Int_t k= 0; k < fNdim; k++) g(k) = random.Gaus(0.,1.);
478
479 // Multiply this vector by Lt to introduce the appropriate correlations
480 g *= (*Lt);
481
482 // Add the mean value offsets and store the results
483 for (int j=1; j<=fNdim; j++) {
484 fToyhisto->SetBinContent( j, fBdat->GetBinContent(j)+g(j-1) );
485 fToyhisto->SetBinError ( j, fBdat->GetBinError(j) );
486 }
487 unfres = Unfold(GetKReg());
488
489 for (Int_t j=1; j<=fNdim; j++) {
490 for (Int_t k=1; k<=fNdim; k++) {
491 unfcov->SetBinContent(j,k,unfcov->GetBinContent(j,k) + ( (unfres->GetBinContent(j) - toymean->GetBinContent(j))* (unfres->GetBinContent(k) - toymean->GetBinContent(k))/(ntoys-1)) );
492 }
493 }
494 delete unfres;
495 unfres = nullptr;
496 }
497 delete Lt;
498 delete toymean;
500
501 return unfcov;
502}
503
504////////////////////////////////////////////////////////////////////////////////
505/// Determine covariance matrix of unfolded spectrum from finite statistics in
506/// response matrix using pseudo experiments
507/// "ntoys" - number of pseudo experiments used for the propagation
508/// "seed" - seed for pseudo experiments
509
511{
512 fMatToyMode = true;
513 TH1D* unfres = nullptr;
514 TH2D* unfcov = (TH2D*)fAdet->Clone("unfcovmat");
515 unfcov->SetTitle("Toy covariance matrix");
516 for(int i=1; i<=fNdim; i++)
517 for(int j=1; j<=fNdim; j++)
518 unfcov->SetBinContent(i,j,0.);
519
520 //Now the toys for the detector response matrix
521 TRandom3 random(seed);
522
523 fToymat = (TH2D*)fAdet->Clone("toymat");
524 TH1D *toymean = (TH1D*)fXini->Clone("toymean");
525 for (Int_t j=1; j<=fNdim; j++) toymean->SetBinContent(j,0.);
526
527 for (int i=1; i<=ntoys; i++) {
528 for (Int_t k=1; k<=fNdim; k++) {
529 for (Int_t m=1; m<=fNdim; m++) {
530 if (fAdet->GetBinContent(k,m)) {
531 fToymat->SetBinContent(k, m, random.Poisson(fAdet->GetBinContent(k,m)));
532 }
533 }
534 }
535
536 unfres = Unfold(GetKReg());
537
538 for (Int_t j=1; j<=fNdim; j++) {
539 toymean->SetBinContent(j, toymean->GetBinContent(j) + unfres->GetBinContent(j)/ntoys);
540 }
541 delete unfres;
542 unfres = nullptr;
543 }
544
545 // Reset the random seed
546 random.SetSeed(seed);
547
548 for (int i=1; i<=ntoys; i++) {
549 for (Int_t k=1; k<=fNdim; k++) {
550 for (Int_t m=1; m<=fNdim; m++) {
551 if (fAdet->GetBinContent(k,m))
552 fToymat->SetBinContent(k, m, random.Poisson(fAdet->GetBinContent(k,m)));
553 }
554 }
555
556 unfres = Unfold(GetKReg());
557
558 for (Int_t j=1; j<=fNdim; j++) {
559 for (Int_t k=1; k<=fNdim; k++) {
560 unfcov->SetBinContent(j,k,unfcov->GetBinContent(j,k) + ( (unfres->GetBinContent(j) - toymean->GetBinContent(j))*(unfres->GetBinContent(k) - toymean->GetBinContent(k))/(ntoys-1)) );
561 }
562 }
563 delete unfres;
564 unfres = nullptr;
565 }
566 delete toymean;
568
569 return unfcov;
570}
571
572////////////////////////////////////////////////////////////////////////////////
573/// Returns d vector (for choosing appropriate regularisation)
574
576{
577 for (int i=1; i<=fDHist->GetNbinsX(); i++) {
579 }
580 return fDHist;
581}
582
583////////////////////////////////////////////////////////////////////////////////
584/// Returns singular values vector
585
587{
588 return fSVHist;
589}
590
591////////////////////////////////////////////////////////////////////////////////
592/// Returns the computed regularized covariance matrix corresponding to total uncertainties on measured spectrum as passed in the constructor.
593/// Note that this covariance matrix will not contain the effects of forced normalization if spectrum is normalized to unit area.
594
596{
597 return fXtau;
598}
599
600////////////////////////////////////////////////////////////////////////////////
601/// Returns the computed inverse of the covariance matrix
602
604{
605 return fXinv;
606}
607
608////////////////////////////////////////////////////////////////////////////////
609/// Returns the covariance matrix
610
612{
613 return fBcov;
614}
615
616////////////////////////////////////////////////////////////////////////////////
617/// Fill 1D histogram into vector
618
619void TSVDUnfold::H2V( const TH1D* histo, TVectorD& vec )
620{
621 for (Int_t i=0; i<histo->GetNbinsX(); i++) vec(i) = histo->GetBinContent(i+1);
622}
623
624////////////////////////////////////////////////////////////////////////////////
625/// Fill 1D histogram errors into vector
626
627void TSVDUnfold::H2Verr( const TH1D* histo, TVectorD& vec )
628{
629 for (Int_t i=0; i<histo->GetNbinsX(); i++) vec(i) = histo->GetBinError(i+1);
630}
631
632////////////////////////////////////////////////////////////////////////////////
633/// Fill vector into 1D histogram
634
635void TSVDUnfold::V2H( const TVectorD& vec, TH1D& histo )
636{
637 for(Int_t i=0; i<vec.GetNrows(); i++) histo.SetBinContent(i+1, vec(i));
638}
639
640////////////////////////////////////////////////////////////////////////////////
641/// Fill 2D histogram into matrix
642
643void TSVDUnfold::H2M( const TH2D* histo, TMatrixD& mat )
644{
645 for (Int_t j=0; j<histo->GetNbinsX(); j++) {
646 for (Int_t i=0; i<histo->GetNbinsY(); i++) {
647 mat(i,j) = histo->GetBinContent(i+1,j+1);
648 }
649 }
650}
651
652////////////////////////////////////////////////////////////////////////////////
653/// Fill 2D histogram into matrix
654
655void TSVDUnfold::M2H( const TMatrixD& mat, TH2D& histo )
656{
657 for (Int_t j=0; j<mat.GetNcols(); j++) {
658 for (Int_t i=0; i<mat.GetNrows(); i++) {
659 histo.SetBinContent(i+1,j+1, mat(i,j));
660 }
661 }
662}
663
664////////////////////////////////////////////////////////////////////////////////
665/// Divide entries of two vectors
666
668{
669 TVectorD quot(vec1.GetNrows());
670 for (Int_t i=0; i<vec1.GetNrows(); i++) {
671 if (vec2(i) != 0) quot(i) = vec1(i) / vec2(i);
672 else {
673 if (zero) quot(i) = 0;
674 else quot(i) = vec1(i);
675 }
676 }
677 return quot;
678}
679
680////////////////////////////////////////////////////////////////////////////////
681/// Divide matrix entries by vector
682
684{
685 TMatrixD quotmat(mat.GetNrows(), mat.GetNcols());
686 for (Int_t i=0; i<mat.GetNrows(); i++) {
687 for (Int_t j=0; j<mat.GetNcols(); j++) {
688 if (vec(i) != 0) quotmat(i,j) = mat(i,j) / vec(i);
689 else {
690 if (zero) quotmat(i,j) = 0;
691 else quotmat(i,j) = mat(i,j);
692 }
693 }
694 }
695 return quotmat;
696}
697
698////////////////////////////////////////////////////////////////////////////////
699/// Multiply entries of two vectors
700
702{
703 TVectorD res(vec1.GetNrows());
704 for (Int_t i=0; i<vec1.GetNrows(); i++) res(i) = vec1(i) * vec2(i);
705 return res;
706}
707
708////////////////////////////////////////////////////////////////////////////////
709/// Compute curvature of vector
710
712{
713 return vec*(curv*vec);
714}
715
716////////////////////////////////////////////////////////////////////////////////
717
719{
720 Double_t eps = 0.00001;
721
722 Int_t ndim = tCurv.GetNrows();
723
724 // Init
725 tCurv *= 0;
726 tC *= 0;
727
728 if (fDdim == 0) for (Int_t i=0; i<ndim; i++) tC(i,i) = 1;
729 else if (fDdim == 1) {
730 for (Int_t i=0; i<ndim; i++) {
731 if (i < ndim-1) tC(i,i+1) = 1.0;
732 tC(i,i) = 1.0;
733 }
734 }
735 else if (fDdim == 2) {
736 for (Int_t i=0; i<ndim; i++) {
737 if (i > 0) tC(i,i-1) = 1.0;
738 if (i < ndim-1) tC(i,i+1) = 1.0;
739 tC(i,i) = -2.0;
740 }
741 tC(0,0) = -1.0;
742 tC(ndim-1,ndim-1) = -1.0;
743 }
744 else if (fDdim == 3) {
745 for (Int_t i=1; i<ndim-2; i++) {
746 tC(i,i-1) = 1.0;
747 tC(i,i) = -3.0;
748 tC(i,i+1) = 3.0;
749 tC(i,i+2) = -1.0;
750 }
751 }
752 else if (fDdim==4) {
753 for (Int_t i=0; i<ndim; i++) {
754 if (i > 0) tC(i,i-1) = -4.0;
755 if (i < ndim-1) tC(i,i+1) = -4.0;
756 if (i > 1) tC(i,i-2) = 1.0;
757 if (i < ndim-2) tC(i,i+2) = 1.0;
758 tC(i,i) = 6.0;
759 }
760 tC(0,0) = 2.0;
761 tC(ndim-1,ndim-1) = 2.0;
762 tC(0,1) = -3.0;
763 tC(ndim-2,ndim-1) = -3.0;
764 tC(1,0) = -3.0;
765 tC(ndim-1,ndim-2) = -3.0;
766 tC(1,1) = 6.0;
767 tC(ndim-2,ndim-2) = 6.0;
768 }
769 else if (fDdim == 5) {
770 for (Int_t i=2; i < ndim-3; i++) {
771 tC(i,i-2) = 1.0;
772 tC(i,i-1) = -5.0;
773 tC(i,i) = 10.0;
774 tC(i,i+1) = -10.0;
775 tC(i,i+2) = 5.0;
776 tC(i,i+3) = -1.0;
777 }
778 }
779 else if (fDdim == 6) {
780 for (Int_t i = 3; i < ndim - 3; i++) {
781 tC(i,i-3) = 1.0;
782 tC(i,i-2) = -6.0;
783 tC(i,i-1) = 15.0;
784 tC(i,i) = -20.0;
785 tC(i,i+1) = 15.0;
786 tC(i,i+2) = -6.0;
787 tC(i,i+3) = 1.0;
788 }
789 }
790
791 // Add epsilon to avoid singularities
792 for (Int_t i=0; i<ndim; i++) tC(i,i) = tC(i,i) + eps;
793
794 //Get curvature matrix
795 for (Int_t i=0; i<ndim; i++) {
796 for (Int_t j=0; j<ndim; j++) {
797 for (Int_t k=0; k<ndim; k++) {
798 tCurv(i,j) = tCurv(i,j) + tC(k,i)*tC(k,j);
799 }
800 }
801 }
802}
803
804////////////////////////////////////////////////////////////////////////////////
805
807{
808 fDHist = new TH1D( "dd", "d vector after orthogonal transformation", fNdim, 0, fNdim );
809 fDHist->Sumw2();
810
811 fSVHist = new TH1D( "sv", "Singular values of AC^-1", fNdim, 0, fNdim );
812 fSVHist->Sumw2();
813
814 fXtau = (TH2D*)fAdet->Clone("Xtau");
815 fXtau->SetTitle("Regularized covariance matrix");
816 fXtau->Sumw2();
817
818 fXinv = (TH2D*)fAdet->Clone("Xinv");
819 fXinv->SetTitle("Inverse covariance matrix");
820 fXinv->Sumw2();
821}
822
823////////////////////////////////////////////////////////////////////////////////
824/// naive regularised inversion cuts off small elements
825
827{
828 // init reduced matrix
829 const UInt_t n = mat.GetNrows();
830 UInt_t nn = 0;
831
832 UInt_t *ipos = new UInt_t[n];
833 // UInt_t ipos[n];
834
835 // find max diagonal entries
836 Double_t ymax = 0;
837 for (UInt_t i=0; i<n; i++) if (TMath::Abs(mat[i][i]) > ymax) ymax = TMath::Abs(mat[i][i]);
838
839 for (UInt_t i=0; i<n; i++) {
840
841 // save position of accepted entries
842 if (TMath::Abs(mat[i][i])/ymax > eps) ipos[nn++] = i;
843 }
844
845 // effective matrix
847 for (UInt_t in=0; in<nn; in++) for (UInt_t jn=0; jn<nn; jn++) matwork(in,jn) = 0;
848
849 // fill non-zero effective working matrix
850 for (UInt_t in=0; in<nn; in++) {
851
852 matwork[in][in] = mat[ipos[in]][ipos[in]];
853 for (UInt_t jn=in+1; jn<nn; jn++) {
854 matwork[in][jn] = mat[ipos[in]][ipos[jn]];
855 matwork[jn][in] = matwork[in][jn];
856 }
857 }
858
859 // invert
860 matwork.Invert();
861
862 // reinitialise old matrix
863 for (UInt_t i=0; i<n; i++) for (UInt_t j=0; j<n; j++) mat[i][j] = 0;
864
865 // refill inverted matrix in old one
866 for (UInt_t in=0; in<nn; in++) {
867 mat[ipos[in]][ipos[in]] = matwork[in][in];
868 for (UInt_t jn=in+1; jn<nn; jn++) {
869 mat[ipos[in]][ipos[jn]] = matwork[in][jn];
870 mat[ipos[jn]][ipos[in]] = mat[ipos[in]][ipos[jn]];
871 }
872 }
873 delete [] ipos;
874}
875
876////////////////////////////////////////////////////////////////////////////////
877/// Helper routine to compute chi-squared between distributions using the computed inverse of the covariance matrix for the unfolded spectrum as given in paper.
878
880{
881 UInt_t n = truspec.GetNbinsX();
882
883 // compute chi2
884 Double_t chi2 = 0;
885 for (UInt_t i=0; i<n; i++) {
886 for (UInt_t j=0; j<n; j++) {
887 chi2 += ( (truspec.GetBinContent( i+1 )-unfspec.GetBinContent( i+1 )) *
888 (truspec.GetBinContent( j+1 )-unfspec.GetBinContent( j+1 )) * fXinv->GetBinContent(i+1,j+1) );
889 }
890 }
891
892 return chi2;
893}
894
#define g(i)
Definition RSha256.hxx:105
#define a(i)
Definition RSha256.hxx:99
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
float ymax
TMatrixT< Double_t > TMatrixD
Definition TMatrixDfwd.h:23
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2495
Single Value Decomposition class.
Definition TDecompSVD.h:24
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:927
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6766
virtual Int_t GetNbinsY() const
Definition TH1.h:543
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition TH1.cxx:9101
virtual Int_t GetNbinsX() const
Definition TH1.h:542
virtual void SetBinError(Int_t bin, Double_t error)
Set the bin Error Note that this resets the bin eror option to be of Normal Type and for the non-empt...
Definition TH1.cxx:9244
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition TH1.cxx:9260
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5076
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition TH1.cxx:9061
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:356
void SetBinContent(Int_t bin, Double_t content) override
Set bin content.
Definition TH2.cxx:2578
Double_t GetBinContent(Int_t binx, Int_t biny) const override
Definition TH2.h:97
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1099
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:1045
Random number generator class based on M.
Definition TRandom3.h:27
SVD Approach to Data Unfolding.
Definition TSVDUnfold.h:46
TH1D * GetSV() const
Returns singular values vector.
TH2D * GetBCov() const
Returns the covariance matrix.
TH1D * fSVHist
! Distribution of singular values
Definition TSVDUnfold.h:133
TH2D * GetXtau() const
Returns the computed regularized covariance matrix corresponding to total uncertainties on measured s...
Bool_t fToyMode
! Internal switch for covariance matrix propagation
Definition TSVDUnfold.h:151
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.
static void RegularisedSymMatInvert(TMatrixDSym &mat, Double_t eps=1e-3)
naive regularised inversion cuts off small elements
static TVectorD CompProd(const TVectorD &vec1, const TVectorD &vec2)
Multiply entries of two vectors.
static void M2H(const TMatrixD &mat, TH2D &histo)
Fill 2D histogram into matrix.
Bool_t fMatToyMode
! Internal switch for evaluation of statistical uncertainties from response matrix
Definition TSVDUnfold.h:152
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,...
Bool_t fNormalize
! Normalize unfolded spectrum to 1
Definition TSVDUnfold.h:130
const TH1D * fBini
Reconstructed distribution (MC)
Definition TSVDUnfold.h:142
TH1D * Unfold(Int_t kreg)
Perform the unfolding with regularisation parameter kreg.
Double_t ComputeChiSquared(const TH1D &truspec, const TH1D &unfspec)
Helper routine to compute chi-squared between distributions using the computed inverse of the covaria...
TH2D * fToymat
! Toy MC detector response matrix
Definition TSVDUnfold.h:150
static Double_t GetCurvature(const TVectorD &vec, const TMatrixD &curv)
Compute curvature of vector.
const TH1D * fBdat
Measured distribution (data)
Definition TSVDUnfold.h:140
~TSVDUnfold() override
Destructor.
static TMatrixD MatDivVec(const TMatrixD &mat, const TVectorD &vec, Int_t zero=0)
Divide matrix entries by vector.
TH2D * fXinv
! Computed inverse of covariance matrix
Definition TSVDUnfold.h:135
static TVectorD VecDiv(const TVectorD &vec1, const TVectorD &vec2, Int_t zero=0)
Divide entries of two vectors.
void FillCurvatureMatrix(TMatrixD &tCurv, TMatrixD &tC) const
TH1D * GetD() const
Returns d vector (for choosing appropriate regularisation)
Int_t GetKReg() const
Definition TSVDUnfold.h:86
TH2D * fXtau
! Computed regularized covariance matrix
Definition TSVDUnfold.h:134
TH1D * fDHist
! Distribution of d (for checking regularization)
Definition TSVDUnfold.h:132
static void H2Verr(const TH1D *histo, TVectorD &vec)
Fill 1D histogram errors into vector.
Int_t fDdim
! Derivative for curvature matrix
Definition TSVDUnfold.h:129
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 * GetAdetCovMatrix(Int_t ntoys, Int_t seed=1)
Determine covariance matrix of unfolded spectrum from finite statistics in response matrix using pseu...
const TH2D * fAdet
Detector response matrix.
Definition TSVDUnfold.h:144
TH2D * GetXinv() const
Returns the computed inverse of the covariance matrix.
static void H2V(const TH1D *histo, TVectorD &vec)
Fill 1D histogram into vector.
Int_t fKReg
! Regularisation parameter
Definition TSVDUnfold.h:131
TH1D * fToyhisto
! Toy MC histogram
Definition TSVDUnfold.h:149
void InitHistos()
TH2D * fBcov
Covariance matrix of measured distribution (data)
Definition TSVDUnfold.h:141
const TH1D * fXini
Truth distribution (MC)
Definition TSVDUnfold.h:143
Int_t fNdim
! Truth and reconstructed dimensions
Definition TSVDUnfold.h:128
Basic string class.
Definition TString.h:138
const Int_t n
Definition legend1.C:16
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:732
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124
TMarker m
Definition textangle.C:8