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