Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TLimit.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Christophe.Delaere@cern.ch 21/08/2002
3
4///////////////////////////////////////////////////////////////////////////
5//
6// TLimit
7//
8// Class to compute 95% CL limits
9//
10// adapted from the mclimit code from Tom Junk (CLs method)
11// see http://root.cern/root/doc/TomJunk.pdf
12// see http://cern.ch/thomasj/searchlimits/ecl.html
13// see: Tom Junk,NIM A434, p. 435-443, 1999
14//
15// see also the following interesting references:
16// Alex Read, "Presentation of search results: the CLs technique"
17// Journal of Physics G: Nucl. Part. Phys. 28 2693-2704 (2002).
18// http://www.iop.org/EJ/abstract/0954-3899/28/10/313/
19//
20// A nice article is also available in the CERN yellow report with the proceeding
21// of the 2000 CERN workshop on confidence intervals.
22//
23// Alex Read, "Modified Frequentist Analysis of Search Results (The CLs Method)"
24// CERN 2000-005 (30 May 2000)
25//
26///////////////////////////////////////////////////////////////////////////
27
28/** \class TLimit
29 \legacy{TLimit, Consider switching to RooStats.}
30 \ingroup Hist
31 Algorithm to compute 95% CL limits using the Likelihood ratio semi-bayesian method.
32
33 Implemented by C. Delaere from the mclimit code written by Tom Junk [HEP-EX/9902006].
34 See [http://cern.ch/thomasj/searchlimits/ecl.html](http://cern.ch/thomasj/searchlimits/ecl.html) for more details.
35
36 It takes signal, background and data histograms wrapped in a
37 TLimitDataSource as input and runs a set of Monte Carlo experiments in
38 order to compute the limits. If needed, inputs are fluctuated according
39 to systematics. The output is a TConfidenceLevel.
40
41 The class TLimitDataSource takes the signal, background and data histograms as well as different
42 systematics sources to form the TLimit input.
43
44 The class TConfidenceLevel represents the final result of the TLimit algorithm. It is created just after the
45 time-consuming part and can be stored in a TFile for further processing.
46 It contains light methods to return CLs, CLb and other interesting
47 quantities.
48
49 The actual algorithm...
50
51 From an input (TLimitDataSource) it produces an output TConfidenceLevel.
52 For this, nmc Monte Carlo experiments are performed.
53 As usual, the larger this number, the longer the compute time,
54 but the better the result.
55
56 Supposing that there is a plotfile.root file containing 3 histograms
57 (signal, background and data), you can imagine doing things like:
58
59~~~{.cpp}
60 TFile* infile=new TFile("plotfile.root","READ");
61 infile->cd();
62 TH1* sh=(TH1*)infile->Get("signal");
63 TH1* bh=(TH1*)infile->Get("background");
64 TH1* dh=(TH1*)infile->Get("data");
65 TLimitDataSource* mydatasource = new TLimitDataSource(sh,bh,dh);
66 TConfidenceLevel *myconfidence = TLimit::ComputeLimit(mydatasource,50000);
67 std::cout << " CLs : " << myconfidence->CLs() << std::endl;
68 std::cout << " CLsb : " << myconfidence->CLsb() << std::endl;
69 std::cout << " CLb : " << myconfidence->CLb() << std::endl;
70 std::cout << "< CLs > : " << myconfidence->GetExpectedCLs_b() << std::endl;
71 std::cout << "< CLsb > : " << myconfidence->GetExpectedCLsb_b() << std::endl;
72 std::cout << "< CLb > : " << myconfidence->GetExpectedCLb_b() << std::endl;
73 delete myconfidence;
74 delete mydatasource;
75 infile->Close();
76~~~
77 More information can still be found on [this page](http://cern.ch/aleph-proj-alphapp/doc/tlimit.html)
78 \see https://doi.org/10.1088/0954-3899/28/10/313 and https://cds.cern.ch/record/451614/files/open-2000-205.pdf
79 \note see note about: "Should I use TRolke, TFeldmanCousins, TLimit?" in the TRolke class description.
80 */
81
82#include "TLimit.h"
83#include "TArrayD.h"
84#include "TVectorD.h"
85#include "TOrdCollection.h"
86#include "TConfidenceLevel.h"
87#include "TLimitDataSource.h"
88#include "TRandom3.h"
89#include "TH1.h"
90#include "TObjArray.h"
91#include "TMath.h"
92#include "TIterator.h"
93#include "TObjString.h"
94#include "TClassTable.h"
95#include "Riostream.h"
96
98
101
103 Int_t nmc, bool stat,
104 TRandom * generator)
105{
106 // The final object returned...
108 // The random generator used...
109 TRandom *myrandom = generator ? generator : new TRandom3;
110 // Compute some total quantities on all the channels
111 Int_t maxbins = 0;
112 Double_t nsig = 0;
113 Double_t nbg = 0;
114 Int_t ncand = 0;
115 Int_t i;
116 for (i = 0; i <= data->GetSignal()->GetLast(); i++) {
117 maxbins = (((TH1 *) (data->GetSignal()->At(i)))->GetNbinsX() + 2) > maxbins ?
118 (((TH1 *) (data->GetSignal()->At(i)))->GetNbinsX() + 2) : maxbins;
119 nsig += ((TH1 *) (data->GetSignal()->At(i)))->Integral();
120 nbg += ((TH1 *) (data->GetBackground()->At(i)))->Integral();
121 ncand += (Int_t) ((TH1 *) (data->GetCandidates()->At(i)))->Integral();
122 }
123 result->SetBtot(nbg);
124 result->SetStot(nsig);
125 result->SetDtot(ncand);
126 Double_t buffer = 0;
127 fgTable->Set(maxbins * (data->GetSignal()->GetLast() + 1));
128 for (Int_t channel = 0; channel <= data->GetSignal()->GetLast(); channel++)
129 for (Int_t bin = 0;
130 bin <= ((TH1 *) (data->GetSignal()->At(channel)))->GetNbinsX()+1;
131 bin++) {
132 Double_t s = (Double_t) ((TH1 *) (data->GetSignal()->At(channel)))->GetBinContent(bin);
133 Double_t b = (Double_t) ((TH1 *) (data->GetBackground()->At(channel)))->GetBinContent(bin);
134 Double_t d = (Double_t) ((TH1 *) (data->GetCandidates()->At(channel)))->GetBinContent(bin);
135 // Compute the value of the "-2lnQ" for the actual data
136 if ((b == 0) && (s > 0)) {
137 std::cout << "WARNING: Ignoring bin " << bin << " of channel "
138 << channel << " which has s=" << s << " but b=" << b << std::endl;
139 std::cout << " Maybe the MC statistic has to be improved..." << std::endl;
140 }
141 if ((s > 0) && (b > 0))
142 buffer += LogLikelihood(s, b, b, d);
143 // precompute the log(1+s/b)'s in an array to speed up computation
144 // background-free bins are set to have a maximum t.s. value
145 // for protection (corresponding to s/b of about 5E8)
146 if ((s > 0) && (b > 0))
147 fgTable->AddAt(LogLikelihood(s, b, b, 1), (channel * maxbins) + bin);
148 else if ((s > 0) && (b == 0))
149 fgTable->AddAt(20, (channel * maxbins) + bin);
150 }
151 result->SetTSD(buffer);
152 // accumulate MC experiments. Hold the test statistic function fixed, but
153 // fluctuate s and b within syst. errors for computing probabilities of
154 // having that outcome. (Alex Read's prescription -- errors are on the ensemble,
155 // not on the observed test statistic. This technique does not split outcomes.)
156 // keep the tstats as sum log(1+s/b). convert to -2lnQ when preparing the results
157 // (reason -- like to keep the < signs right)
158 Double_t *tss = new Double_t[nmc];
159 Double_t *tsb = new Double_t[nmc];
160 Double_t *lrs = new Double_t[nmc];
161 Double_t *lrb = new Double_t[nmc];
162 TLimitDataSource* tmp_src = (TLimitDataSource*)(data->Clone());
163 TLimitDataSource* tmp_src2 = (TLimitDataSource*)(data->Clone());
164 for (i = 0; i < nmc; i++) {
165 tss[i] = 0;
166 tsb[i] = 0;
167 lrs[i] = 0;
168 lrb[i] = 0;
169 // fluctuate signal and background
170 TLimitDataSource* fluctuated = Fluctuate(data, tmp_src, !i, myrandom, stat) ? tmp_src : data;
171 // make an independent set of fluctuations for reweighting pseudoexperiments
172 // it turns out using the same fluctuated ones for numerator and denominator
173 // is biased.
174 TLimitDataSource* fluctuated2 = Fluctuate(data, tmp_src2, false, myrandom, stat) ? tmp_src2 : data;
175
176 for (Int_t channel = 0;
177 channel <= fluctuated->GetSignal()->GetLast(); channel++) {
178 for (Int_t bin = 0;
179 bin <=((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetNbinsX()+1;
180 bin++) {
181 if ((Double_t) ((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetBinContent(bin) != 0) {
182 // s+b hypothesis
183 Double_t rate = (Double_t) ((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetBinContent(bin) +
184 (Double_t) ((TH1 *) (fluctuated->GetBackground()->At(channel)))->GetBinContent(bin);
185 Double_t rand = myrandom->Poisson(rate);
186 tss[i] += rand * fgTable->At((channel * maxbins) + bin);
187 Double_t s = (Double_t) ((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetBinContent(bin);
188 Double_t s2= (Double_t) ((TH1 *) (fluctuated2->GetSignal()->At(channel)))->GetBinContent(bin);
189 Double_t b = (Double_t) ((TH1 *) (fluctuated->GetBackground()->At(channel)))->GetBinContent(bin);
190 Double_t b2= (Double_t) ((TH1 *) (fluctuated2->GetBackground()->At(channel)))->GetBinContent(bin);
191 if ((s > 0) && (b2 > 0))
192 lrs[i] += LogLikelihood(s, b, b2, rand) - s - b + b2;
193 else if ((s > 0) && (b2 == 0))
194 lrs[i] += 20 * rand - s;
195 // b hypothesis
196 rate = (Double_t) ((TH1 *) (fluctuated->GetBackground()->At(channel)))->GetBinContent(bin);
197 rand = myrandom->Poisson(rate);
198 tsb[i] += rand * fgTable->At((channel * maxbins) + bin);
199 if ((s2 > 0) && (b > 0))
200 lrb[i] += LogLikelihood(s2, b2, b, rand) - s2 - b2 + b;
201 else if ((s > 0) && (b == 0))
202 lrb[i] += 20 * rand - s;
203 }
204 }
205 }
206 lrs[i] = TMath::Exp(lrs[i] < 710 ? lrs[i] : 710);
207 lrb[i] = TMath::Exp(lrb[i] < 710 ? lrb[i] : 710);
208 }
209 delete tmp_src;
210 delete tmp_src2;
211 // lrs and lrb are the LR's (no logs) = prob(s+b)/prob(b) for
212 // that choice of s and b within syst. errors in the ensemble. These are
213 // the MC experiment weights for relating the s+b and b PDF's of the unsmeared
214 // test statistic (in which cas one can use another test statistic if one likes).
215
216 // Now produce the output object.
217 // The final quantities are computed on-demand form the arrays tss, tsb, lrs and lrb.
218 result->SetTSS(tss);
219 result->SetTSB(tsb);
220 result->SetLRS(lrs);
221 result->SetLRB(lrb);
222 if (!generator)
223 delete myrandom;
224 return result;
225}
226
228 bool init, TRandom * generator, bool stat)
229{
230 // initialisation: create a sorted list of all the names of systematics
231 if (init) {
232 // create a "map" with the systematics names
233 TIter errornames = input->GetErrorNames()->MakeIterator();
234 TObjArray *listofnames = nullptr;
235 delete fgSystNames;
237 while ((listofnames = ((TObjArray *) errornames.Next()))) {
238 TObjString *name = nullptr;
239 TIter loniter = listofnames->MakeIterator();
240 while ((name = (TObjString *) loniter.Next()))
241 if ((fgSystNames->IndexOf(name)) < 0)
243 }
244 fgSystNames->Sort();
245 }
246 // if the output is not given, create it from the input
247 if (!output)
248 output = (TLimitDataSource*)(input->Clone());
249 // if there are no systematics, just returns the input as "fluctuated" output
250 if ((fgSystNames->GetSize() <= 0)&&(!stat)) {
251 return false;
252 }
253 // if there are only stat, just fluctuate stats.
254 if (fgSystNames->GetSize() <= 0) {
255 output->SetOwner();
256 for (Int_t channel = 0; channel <= input->GetSignal()->GetLast(); channel++) {
257 TH1 *newsignal = (TH1*)(output->GetSignal()->At(channel));
258 TH1 *oldsignal = (TH1*)(input->GetSignal()->At(channel));
259 if(stat)
260 for(int i=1; i<=newsignal->GetNbinsX(); i++) {
261 newsignal->SetBinContent(i,oldsignal->GetBinContent(i)+generator->Gaus(0,oldsignal->GetBinError(i)));
262 }
263 newsignal->SetDirectory(nullptr);
264 TH1 *newbackground = (TH1*)(output->GetBackground()->At(channel));
265 TH1 *oldbackground = (TH1*)(input->GetBackground()->At(channel));
266 if(stat)
267 for(int i=1; i<=newbackground->GetNbinsX(); i++)
268 newbackground->SetBinContent(i,oldbackground->GetBinContent(i)+generator->Gaus(0,oldbackground->GetBinError(i)));
269 newbackground->SetDirectory(nullptr);
270 }
271 return true;
272 }
273 // Find a choice for the random variation and
274 // re-toss all random numbers if any background or signal
275 // goes negative. (background = 0 is bad too, so put a little protection
276 // around it -- must have at least 10% of the bg estimate).
277 Bool_t retoss = kTRUE;
278 Double_t *serrf = new Double_t[(input->GetSignal()->GetLast()) + 1];
279 Double_t *berrf = new Double_t[(input->GetSignal()->GetLast()) + 1];
280 do {
281 Double_t *toss = new Double_t[fgSystNames->GetSize()];
282 for (Int_t i = 0; i < fgSystNames->GetSize(); i++)
283 toss[i] = generator->Gaus(0, 1);
284 retoss = kFALSE;
285 for (Int_t channel = 0;
286 channel <= input->GetSignal()->GetLast();
287 channel++) {
288 serrf[channel] = 0;
289 berrf[channel] = 0;
290 for (Int_t bin = 0;
291 bin <((TVectorD *) (input->GetErrorOnSignal()->At(channel)))->GetNrows();
292 bin++) {
293 serrf[channel] += ((TVectorD *) (input->GetErrorOnSignal()->At(channel)))->operator[](bin) *
294 toss[fgSystNames->BinarySearch((TObjString*) (((TObjArray *) (input->GetErrorNames()->At(channel)))->At(bin)))];
295 berrf[channel] += ((TVectorD *) (input->GetErrorOnBackground()->At(channel)))->operator[](bin) *
296 toss[fgSystNames->BinarySearch((TObjString*) (((TObjArray *) (input->GetErrorNames()->At(channel)))->At(bin)))];
297 }
298 if ((serrf[channel] < -1.0) || (berrf[channel] < -0.9)) {
299 retoss = kTRUE;
300 continue;
301 }
302 }
303 delete[]toss;
304 } while (retoss);
305 // adjust the fluctuated signal and background counts with a legal set
306 // of random fluctuations above.
307 output->SetOwner();
308 for (Int_t channel = 0; channel <= input->GetSignal()->GetLast();
309 channel++) {
310 TH1 *newsignal = (TH1*)(output->GetSignal()->At(channel));
311 TH1 *oldsignal = (TH1*)(input->GetSignal()->At(channel));
312 if(stat)
313 for(int i=1; i<=newsignal->GetNbinsX(); i++)
314 newsignal->SetBinContent(i,oldsignal->GetBinContent(i)+generator->Gaus(0,oldsignal->GetBinError(i)));
315 else
316 for(int i=1; i<=newsignal->GetNbinsX(); i++)
317 newsignal->SetBinContent(i,oldsignal->GetBinContent(i));
318 newsignal->Scale(1 + serrf[channel]);
319 newsignal->SetDirectory(nullptr);
320 TH1 *newbackground = (TH1*)(output->GetBackground()->At(channel));
321 TH1 *oldbackground = (TH1*)(input->GetBackground()->At(channel));
322 if(stat)
323 for(int i=1; i<=newbackground->GetNbinsX(); i++)
324 newbackground->SetBinContent(i,oldbackground->GetBinContent(i)+generator->Gaus(0,oldbackground->GetBinError(i)));
325 else
326 for(int i=1; i<=newbackground->GetNbinsX(); i++)
327 newbackground->SetBinContent(i,oldbackground->GetBinContent(i));
328 newbackground->Scale(1 + berrf[channel]);
329 newbackground->SetDirectory(nullptr);
330 }
331 delete[] serrf;
332 delete[] berrf;
333 return true;
334}
335
337 Int_t nmc, bool stat,
338 TRandom * generator)
339{
340 // Compute limit.
341
342 TLimitDataSource* lds = new TLimitDataSource(s,b,d);
343 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
344 delete lds;
345 return out;
346}
347
349 TVectorD* se, TVectorD* be, TObjArray* l,
350 Int_t nmc, bool stat,
351 TRandom * generator)
352{
353 // Compute limit.
354
355 TLimitDataSource* lds = new TLimitDataSource(s,b,d,se,be,l);
356 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
357 delete lds;
358 return out;
359}
360
362 Int_t nmc,
363 bool stat,
364 TRandom * generator)
365{
366 // Compute limit.
367
368 TH1D* sh = new TH1D("__sh","__sh",1,0,2);
369 sh->Fill(1,s);
370 TH1D* bh = new TH1D("__bh","__bh",1,0,2);
371 bh->Fill(1,b);
372 TH1D* dh = new TH1D("__dh","__dh",1,0,2);
373 dh->Fill(1,d);
374 TLimitDataSource* lds = new TLimitDataSource(sh,bh,dh);
375 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
376 delete lds;
377 delete sh;
378 delete bh;
379 delete dh;
380 return out;
381}
382
384 TVectorD* se, TVectorD* be, TObjArray* l,
385 Int_t nmc,
386 bool stat,
387 TRandom * generator)
388{
389 // Compute limit.
390
391 TH1D* sh = new TH1D("__sh","__sh",1,0,2);
392 sh->Fill(1,s);
393 TH1D* bh = new TH1D("__bh","__bh",1,0,2);
394 bh->Fill(1,b);
395 TH1D* dh = new TH1D("__dh","__dh",1,0,2);
396 dh->Fill(1,d);
397 TLimitDataSource* lds = new TLimitDataSource(sh,bh,dh,se,be,l);
398 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
399 delete lds;
400 delete sh;
401 delete bh;
402 delete dh;
403 return out;
404}
405
407{
408 // Compute LogLikelihood (static function)
409
410 return d*TMath::Log((s+b)/b2);
411}
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
#define ClassImp(name)
Definition Rtypes.h:382
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void input
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
char name[80]
Definition TGX11.cxx:110
Array of doubles (64 bits per element).
Definition TArrayD.h:27
Double_t At(Int_t i) const
Definition TArrayD.h:79
void Set(Int_t n) override
Set size of this array to n doubles.
Definition TArrayD.cxx:106
void AddAt(Double_t c, Int_t i)
Set the double c value at position i in the array.
Definition TArrayD.cxx:94
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
Class to compute 95% CL limits.
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:671
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
virtual void SetDirectory(TDirectory *dir)
By default, when a histogram is created, it is added to the list of histogram objects in the current ...
Definition TH1.cxx:8970
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition TH1.cxx:9096
virtual Int_t GetNbinsX() const
Definition TH1.h:298
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3346
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:9255
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5090
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition TH1.cxx:6633
TObject * Next()
This class serves as input for the TLimit::ComputeLimit method.
virtual TObjArray * GetSignal()
virtual TObjArray * GetBackground()
<div class="legacybox"><h2>Legacy Code</h2> TLimit is a legacy interface: there will be no bug fixes ...
Definition TLimit.h:19
static Double_t LogLikelihood(Double_t s, Double_t b, Double_t b2, Double_t d)
Definition TLimit.cxx:406
static TOrdCollection * fgSystNames
Collection of systematics names.
Definition TLimit.h:51
static TArrayD * fgTable
A log table... just to speed up calculation.
Definition TLimit.h:50
static bool Fluctuate(TLimitDataSource *input, TLimitDataSource *output, bool init, TRandom *, bool stat=false)
Definition TLimit.cxx:227
static TConfidenceLevel * ComputeLimit(TLimitDataSource *data, Int_t nmc=50000, bool stat=false, TRandom *generator=nullptr)
Definition TLimit.cxx:102
An array of TObjects.
Definition TObjArray.h:31
TIterator * MakeIterator(Bool_t dir=kIterForward) const override
Returns an array iterator.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
Collectable string class.
Definition TObjString.h:28
Ordered collection.
void AddLast(TObject *obj) override
Add object at the end of the collection.
void Sort()
If objects in collection are sortable (i.e.
Int_t IndexOf(const TObject *obj) const override
Return index of object in collection.
Int_t BinarySearch(TObject *obj)
Find object using a binary search.
Random number generator class based on M.
Definition TRandom3.h:27
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition TRandom.cxx:275
virtual ULong64_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
Definition TRandom.cxx:404
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
Definition TMath.h:713
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:760
TLine l
Definition textangle.C:4
static void output()