100 TArrayD *TLimit::fgTable =
new TArrayD(0);
101 TOrdCollection *TLimit::fgSystNames =
new TOrdCollection();
103 TConfidenceLevel *TLimit::ComputeLimit(TLimitDataSource * data,
104 Int_t nmc,
bool stat,
108 TConfidenceLevel *result =
new TConfidenceLevel(nmc);
110 TRandom *myrandom = generator ? generator :
new TRandom3;
117 for (i = 0; i <= data->GetSignal()->GetLast(); i++) {
118 maxbins = (((TH1 *) (data->GetSignal()->At(i)))->GetNbinsX() + 2) > maxbins ?
119 (((TH1 *) (data->GetSignal()->At(i)))->GetNbinsX() + 2) : maxbins;
120 nsig += ((TH1 *) (data->GetSignal()->At(i)))->Integral();
121 nbg += ((TH1 *) (data->GetBackground()->At(i)))->Integral();
122 ncand += (Int_t) ((TH1 *) (data->GetCandidates()->At(i)))->Integral();
124 result->SetBtot(nbg);
125 result->SetStot(nsig);
126 result->SetDtot(ncand);
128 fgTable->Set(maxbins * (data->GetSignal()->GetLast() + 1));
129 for (Int_t channel = 0; channel <= data->GetSignal()->GetLast(); channel++)
131 bin <= ((TH1 *) (data->GetSignal()->At(channel)))->GetNbinsX()+1;
133 Double_t s = (Double_t) ((TH1 *) (data->GetSignal()->At(channel)))->GetBinContent(bin);
134 Double_t b = (Double_t) ((TH1 *) (data->GetBackground()->At(channel)))->GetBinContent(bin);
135 Double_t d = (Double_t) ((TH1 *) (data->GetCandidates()->At(channel)))->GetBinContent(bin);
137 if ((b == 0) && (s > 0)) {
138 std::cout <<
"WARNING: Ignoring bin " << bin <<
" of channel "
139 << channel <<
" which has s=" << s <<
" but b=" << b << std::endl;
140 std::cout <<
" Maybe the MC statistic has to be improved..." << std::endl;
142 if ((s > 0) && (b > 0))
143 buffer += LogLikelihood(s, b, b, d);
147 if ((s > 0) && (b > 0))
148 fgTable->AddAt(LogLikelihood(s, b, b, 1), (channel * maxbins) + bin);
149 else if ((s > 0) && (b == 0))
150 fgTable->AddAt(20, (channel * maxbins) + bin);
152 result->SetTSD(buffer);
159 Double_t *tss =
new Double_t[nmc];
160 Double_t *tsb =
new Double_t[nmc];
161 Double_t *lrs =
new Double_t[nmc];
162 Double_t *lrb =
new Double_t[nmc];
163 TLimitDataSource* tmp_src = (TLimitDataSource*)(data->Clone());
164 TLimitDataSource* tmp_src2 = (TLimitDataSource*)(data->Clone());
165 for (i = 0; i < nmc; i++) {
171 TLimitDataSource* fluctuated = Fluctuate(data, tmp_src, !i, myrandom, stat) ? tmp_src : data;
175 TLimitDataSource* fluctuated2 = Fluctuate(data, tmp_src2,
false, myrandom, stat) ? tmp_src2 : data;
177 for (Int_t channel = 0;
178 channel <= fluctuated->GetSignal()->GetLast(); channel++) {
180 bin <=((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetNbinsX()+1;
182 if ((Double_t) ((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetBinContent(bin) != 0) {
184 Double_t rate = (Double_t) ((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetBinContent(bin) +
185 (Double_t) ((TH1 *) (fluctuated->GetBackground()->At(channel)))->GetBinContent(bin);
186 Double_t rand = myrandom->Poisson(rate);
187 tss[i] += rand * fgTable->At((channel * maxbins) + bin);
188 Double_t s = (Double_t) ((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetBinContent(bin);
189 Double_t s2= (Double_t) ((TH1 *) (fluctuated2->GetSignal()->At(channel)))->GetBinContent(bin);
190 Double_t b = (Double_t) ((TH1 *) (fluctuated->GetBackground()->At(channel)))->GetBinContent(bin);
191 Double_t b2= (Double_t) ((TH1 *) (fluctuated2->GetBackground()->At(channel)))->GetBinContent(bin);
192 if ((s > 0) && (b2 > 0))
193 lrs[i] += LogLikelihood(s, b, b2, rand) - s - b + b2;
194 else if ((s > 0) && (b2 == 0))
195 lrs[i] += 20 * rand - s;
197 rate = (Double_t) ((TH1 *) (fluctuated->GetBackground()->At(channel)))->GetBinContent(bin);
198 rand = myrandom->Poisson(rate);
199 tsb[i] += rand * fgTable->At((channel * maxbins) + bin);
200 if ((s2 > 0) && (b > 0))
201 lrb[i] += LogLikelihood(s2, b2, b, rand) - s2 - b2 + b;
202 else if ((s > 0) && (b == 0))
203 lrb[i] += 20 * rand - s;
207 lrs[i] = TMath::Exp(lrs[i] < 710 ? lrs[i] : 710);
208 lrb[i] = TMath::Exp(lrb[i] < 710 ? lrb[i] : 710);
228 bool TLimit::Fluctuate(TLimitDataSource * input, TLimitDataSource * output,
229 bool init, TRandom * generator,
bool stat)
234 TIterator *errornames = input->GetErrorNames()->MakeIterator();
235 TObjArray *listofnames = 0;
237 fgSystNames =
new TOrdCollection();
238 while ((listofnames = ((TObjArray *) errornames->Next()))) {
239 TObjString *name = 0;
240 TIterator *loniter = listofnames->MakeIterator();
241 while ((name = (TObjString *) (loniter->Next())))
242 if ((fgSystNames->IndexOf(name)) < 0)
243 fgSystNames->AddLast(name);
249 output = (TLimitDataSource*)(input->Clone());
251 if ((fgSystNames->GetSize() <= 0)&&(!stat)) {
255 if (fgSystNames->GetSize() <= 0) {
257 for (Int_t channel = 0; channel <= input->GetSignal()->GetLast(); channel++) {
258 TH1 *newsignal = (TH1*)(output->GetSignal()->At(channel));
259 TH1 *oldsignal = (TH1*)(input->GetSignal()->At(channel));
261 for(
int i=1; i<=newsignal->GetNbinsX(); i++) {
262 newsignal->SetBinContent(i,oldsignal->GetBinContent(i)+generator->Gaus(0,oldsignal->GetBinError(i)));
264 newsignal->SetDirectory(0);
265 TH1 *newbackground = (TH1*)(output->GetBackground()->At(channel));
266 TH1 *oldbackground = (TH1*)(input->GetBackground()->At(channel));
268 for(
int i=1; i<=newbackground->GetNbinsX(); i++)
269 newbackground->SetBinContent(i,oldbackground->GetBinContent(i)+generator->Gaus(0,oldbackground->GetBinError(i)));
270 newbackground->SetDirectory(0);
278 Bool_t retoss = kTRUE;
279 Double_t *serrf =
new Double_t[(input->GetSignal()->GetLast()) + 1];
280 Double_t *berrf =
new Double_t[(input->GetSignal()->GetLast()) + 1];
282 Double_t *toss =
new Double_t[fgSystNames->GetSize()];
283 for (Int_t i = 0; i < fgSystNames->GetSize(); i++)
284 toss[i] = generator->Gaus(0, 1);
286 for (Int_t channel = 0;
287 channel <= input->GetSignal()->GetLast();
292 bin <((TVectorD *) (input->GetErrorOnSignal()->At(channel)))->GetNrows();
294 serrf[channel] += ((TVectorD *) (input->GetErrorOnSignal()->At(channel)))->operator[](bin) *
295 toss[fgSystNames->BinarySearch((TObjString*) (((TObjArray *) (input->GetErrorNames()->At(channel)))->At(bin)))];
296 berrf[channel] += ((TVectorD *) (input->GetErrorOnBackground()->At(channel)))->operator[](bin) *
297 toss[fgSystNames->BinarySearch((TObjString*) (((TObjArray *) (input->GetErrorNames()->At(channel)))->At(bin)))];
299 if ((serrf[channel] < -1.0) || (berrf[channel] < -0.9)) {
309 for (Int_t channel = 0; channel <= input->GetSignal()->GetLast();
311 TH1 *newsignal = (TH1*)(output->GetSignal()->At(channel));
312 TH1 *oldsignal = (TH1*)(input->GetSignal()->At(channel));
314 for(
int i=1; i<=newsignal->GetNbinsX(); i++)
315 newsignal->SetBinContent(i,oldsignal->GetBinContent(i)+generator->Gaus(0,oldsignal->GetBinError(i)));
317 for(
int i=1; i<=newsignal->GetNbinsX(); i++)
318 newsignal->SetBinContent(i,oldsignal->GetBinContent(i));
319 newsignal->Scale(1 + serrf[channel]);
320 newsignal->SetDirectory(0);
321 TH1 *newbackground = (TH1*)(output->GetBackground()->At(channel));
322 TH1 *oldbackground = (TH1*)(input->GetBackground()->At(channel));
324 for(
int i=1; i<=newbackground->GetNbinsX(); i++)
325 newbackground->SetBinContent(i,oldbackground->GetBinContent(i)+generator->Gaus(0,oldbackground->GetBinError(i)));
327 for(
int i=1; i<=newbackground->GetNbinsX(); i++)
328 newbackground->SetBinContent(i,oldbackground->GetBinContent(i));
329 newbackground->Scale(1 + berrf[channel]);
330 newbackground->SetDirectory(0);
337 TConfidenceLevel *TLimit::ComputeLimit(TH1* s, TH1* b, TH1* d,
338 Int_t nmc,
bool stat,
343 TLimitDataSource* lds =
new TLimitDataSource(s,b,d);
344 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
349 TConfidenceLevel *TLimit::ComputeLimit(TH1* s, TH1* b, TH1* d,
350 TVectorD* se, TVectorD* be, TObjArray* l,
351 Int_t nmc,
bool stat,
356 TLimitDataSource* lds =
new TLimitDataSource(s,b,d,se,be,l);
357 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
362 TConfidenceLevel *TLimit::ComputeLimit(Double_t s, Double_t b, Int_t d,
369 TH1D* sh =
new TH1D(
"__sh",
"__sh",1,0,2);
371 TH1D* bh =
new TH1D(
"__bh",
"__bh",1,0,2);
373 TH1D* dh =
new TH1D(
"__dh",
"__dh",1,0,2);
375 TLimitDataSource* lds =
new TLimitDataSource(sh,bh,dh);
376 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
384 TConfidenceLevel *TLimit::ComputeLimit(Double_t s, Double_t b, Int_t d,
385 TVectorD* se, TVectorD* be, TObjArray* l,
392 TH1D* sh =
new TH1D(
"__sh",
"__sh",1,0,2);
394 TH1D* bh =
new TH1D(
"__bh",
"__bh",1,0,2);
396 TH1D* dh =
new TH1D(
"__dh",
"__dh",1,0,2);
398 TLimitDataSource* lds =
new TLimitDataSource(sh,bh,dh,se,be,l);
399 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
407 Double_t TLimit::LogLikelihood(Double_t s, Double_t b, Double_t b2, Double_t d)
411 return d*TMath::Log((s+b)/b2);