00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031 #include "TLimit.h"
00032 #include "TArrayD.h"
00033 #include "TVectorD.h"
00034 #include "TOrdCollection.h"
00035 #include "TConfidenceLevel.h"
00036 #include "TLimitDataSource.h"
00037 #include "TRandom3.h"
00038 #include "TH1.h"
00039 #include "TObjArray.h"
00040 #include "TMath.h"
00041 #include "TIterator.h"
00042 #include "TObjString.h"
00043 #include "TClassTable.h"
00044 #include "Riostream.h"
00045
00046 ClassImp(TLimit)
00047
00048 TArrayD *TLimit::fgTable = new TArrayD(0);
00049 TOrdCollection *TLimit::fgSystNames = new TOrdCollection();
00050
00051 TConfidenceLevel *TLimit::ComputeLimit(TLimitDataSource * data,
00052 Int_t nmc, bool stat,
00053 TRandom * generator)
00054 {
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115 TConfidenceLevel *result = new TConfidenceLevel(nmc);
00116
00117 TRandom *myrandom = generator ? generator : new TRandom3;
00118
00119 Int_t maxbins = 0;
00120 Double_t nsig = 0;
00121 Double_t nbg = 0;
00122 Int_t ncand = 0;
00123 Int_t i;
00124 for (i = 0; i <= data->GetSignal()->GetLast(); i++) {
00125 maxbins = (((TH1 *) (data->GetSignal()->At(i)))->GetNbinsX() + 2) > maxbins ?
00126 (((TH1 *) (data->GetSignal()->At(i)))->GetNbinsX() + 2) : maxbins;
00127 nsig += ((TH1 *) (data->GetSignal()->At(i)))->Integral();
00128 nbg += ((TH1 *) (data->GetBackground()->At(i)))->Integral();
00129 ncand += (Int_t) ((TH1 *) (data->GetCandidates()->At(i)))->Integral();
00130 }
00131 result->SetBtot(nbg);
00132 result->SetStot(nsig);
00133 result->SetDtot(ncand);
00134 Double_t buffer = 0;
00135 fgTable->Set(maxbins * (data->GetSignal()->GetLast() + 1));
00136 for (Int_t channel = 0; channel <= data->GetSignal()->GetLast(); channel++)
00137 for (Int_t bin = 0;
00138 bin <= ((TH1 *) (data->GetSignal()->At(channel)))->GetNbinsX()+1;
00139 bin++) {
00140 Double_t s = (Double_t) ((TH1 *) (data->GetSignal()->At(channel)))->GetBinContent(bin);
00141 Double_t b = (Double_t) ((TH1 *) (data->GetBackground()->At(channel)))->GetBinContent(bin);
00142 Double_t d = (Double_t) ((TH1 *) (data->GetCandidates()->At(channel)))->GetBinContent(bin);
00143
00144 if ((b == 0) && (s > 0)) {
00145 cout << "WARNING: Ignoring bin " << bin << " of channel "
00146 << channel << " which has s=" << s << " but b=" << b << endl;
00147 cout << " Maybe the MC statistic has to be improved..." << endl;
00148 }
00149 if ((s > 0) && (b > 0))
00150 buffer += LogLikelihood(s, b, b, d);
00151
00152
00153
00154 if ((s > 0) && (b > 0))
00155 fgTable->AddAt(LogLikelihood(s, b, b, 1), (channel * maxbins) + bin);
00156 else if ((s > 0) && (b == 0))
00157 fgTable->AddAt(20, (channel * maxbins) + bin);
00158 }
00159 result->SetTSD(buffer);
00160
00161
00162
00163
00164
00165
00166 Double_t *tss = new Double_t[nmc];
00167 Double_t *tsb = new Double_t[nmc];
00168 Double_t *lrs = new Double_t[nmc];
00169 Double_t *lrb = new Double_t[nmc];
00170 TLimitDataSource* tmp_src = (TLimitDataSource*)(data->Clone());
00171 TLimitDataSource* tmp_src2 = (TLimitDataSource*)(data->Clone());
00172 for (i = 0; i < nmc; i++) {
00173 tss[i] = 0;
00174 tsb[i] = 0;
00175 lrs[i] = 0;
00176 lrb[i] = 0;
00177
00178 TLimitDataSource* fluctuated = Fluctuate(data, tmp_src, !i, myrandom, stat) ? tmp_src : data;
00179
00180
00181
00182 TLimitDataSource* fluctuated2 = Fluctuate(data, tmp_src2, false, myrandom, stat) ? tmp_src2 : data;
00183
00184 for (Int_t channel = 0;
00185 channel <= fluctuated->GetSignal()->GetLast(); channel++) {
00186 for (Int_t bin = 0;
00187 bin <=((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetNbinsX()+1;
00188 bin++) {
00189 if ((Double_t) ((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetBinContent(bin) != 0) {
00190
00191 Double_t rate = (Double_t) ((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetBinContent(bin) +
00192 (Double_t) ((TH1 *) (fluctuated->GetBackground()->At(channel)))->GetBinContent(bin);
00193 Double_t rand = myrandom->Poisson(rate);
00194 tss[i] += rand * fgTable->At((channel * maxbins) + bin);
00195 Double_t s = (Double_t) ((TH1 *) (fluctuated->GetSignal()->At(channel)))->GetBinContent(bin);
00196 Double_t s2= (Double_t) ((TH1 *) (fluctuated2->GetSignal()->At(channel)))->GetBinContent(bin);
00197 Double_t b = (Double_t) ((TH1 *) (fluctuated->GetBackground()->At(channel)))->GetBinContent(bin);
00198 Double_t b2= (Double_t) ((TH1 *) (fluctuated2->GetBackground()->At(channel)))->GetBinContent(bin);
00199 if ((s > 0) && (b2 > 0))
00200 lrs[i] += LogLikelihood(s, b, b2, rand) - s - b + b2;
00201 else if ((s > 0) && (b2 == 0))
00202 lrs[i] += 20 * rand - s;
00203
00204 rate = (Double_t) ((TH1 *) (fluctuated->GetBackground()->At(channel)))->GetBinContent(bin);
00205 rand = myrandom->Poisson(rate);
00206 tsb[i] += rand * fgTable->At((channel * maxbins) + bin);
00207 if ((s2 > 0) && (b > 0))
00208 lrb[i] += LogLikelihood(s2, b2, b, rand) - s2 - b2 + b;
00209 else if ((s > 0) && (b == 0))
00210 lrb[i] += 20 * rand - s;
00211 }
00212 }
00213 }
00214 lrs[i] = TMath::Exp(lrs[i] < 710 ? lrs[i] : 710);
00215 lrb[i] = TMath::Exp(lrb[i] < 710 ? lrb[i] : 710);
00216 }
00217 delete tmp_src;
00218 delete tmp_src2;
00219
00220
00221
00222
00223
00224
00225
00226 result->SetTSS(tss);
00227 result->SetTSB(tsb);
00228 result->SetLRS(lrs);
00229 result->SetLRB(lrb);
00230 if (!generator)
00231 delete myrandom;
00232 return result;
00233 }
00234
00235 bool TLimit::Fluctuate(TLimitDataSource * input, TLimitDataSource * output,
00236 bool init, TRandom * generator, bool stat)
00237 {
00238
00239 if (init) {
00240
00241 TIterator *errornames = input->GetErrorNames()->MakeIterator();
00242 TObjArray *listofnames = 0;
00243 delete fgSystNames;
00244 fgSystNames = new TOrdCollection();
00245 while ((listofnames = ((TObjArray *) errornames->Next()))) {
00246 TObjString *name = 0;
00247 TIterator *loniter = listofnames->MakeIterator();
00248 while ((name = (TObjString *) (loniter->Next())))
00249 if ((fgSystNames->IndexOf(name)) < 0)
00250 fgSystNames->AddLast(name);
00251 }
00252 fgSystNames->Sort();
00253 }
00254
00255 if (!output)
00256 output = (TLimitDataSource*)(input->Clone());
00257
00258 if ((fgSystNames->GetSize() <= 0)&&(!stat)) {
00259 return 0;
00260 }
00261
00262 if (fgSystNames->GetSize() <= 0) {
00263 output->SetOwner();
00264 for (Int_t channel = 0; channel <= input->GetSignal()->GetLast(); channel++) {
00265 TH1 *newsignal = (TH1*)(output->GetSignal()->At(channel));
00266 TH1 *oldsignal = (TH1*)(input->GetSignal()->At(channel));
00267 if(stat)
00268 for(int i=1; i<=newsignal->GetNbinsX(); i++) {
00269 newsignal->SetBinContent(i,oldsignal->GetBinContent(i)+generator->Gaus(0,oldsignal->GetBinError(i)));
00270 }
00271 newsignal->SetDirectory(0);
00272 TH1 *newbackground = (TH1*)(output->GetBackground()->At(channel));
00273 TH1 *oldbackground = (TH1*)(input->GetBackground()->At(channel));
00274 if(stat)
00275 for(int i=1; i<=newbackground->GetNbinsX(); i++)
00276 newbackground->SetBinContent(i,oldbackground->GetBinContent(i)+generator->Gaus(0,oldbackground->GetBinError(i)));
00277 newbackground->SetDirectory(0);
00278 }
00279 return 1;
00280 }
00281
00282
00283
00284
00285 Bool_t retoss = kTRUE;
00286 Double_t *serrf = new Double_t[(input->GetSignal()->GetLast()) + 1];
00287 Double_t *berrf = new Double_t[(input->GetSignal()->GetLast()) + 1];
00288 do {
00289 Double_t *toss = new Double_t[fgSystNames->GetSize()];
00290 for (Int_t i = 0; i < fgSystNames->GetSize(); i++)
00291 toss[i] = generator->Gaus(0, 1);
00292 retoss = kFALSE;
00293 for (Int_t channel = 0;
00294 channel <= input->GetSignal()->GetLast();
00295 channel++) {
00296 serrf[channel] = 0;
00297 berrf[channel] = 0;
00298 for (Int_t bin = 0;
00299 bin <((TVectorD *) (input->GetErrorOnSignal()->At(channel)))->GetNrows();
00300 bin++) {
00301 serrf[channel] += ((TVectorD *) (input->GetErrorOnSignal()->At(channel)))->operator[](bin) *
00302 toss[fgSystNames->BinarySearch((TObjString*) (((TObjArray *) (input->GetErrorNames()->At(channel)))->At(bin)))];
00303 berrf[channel] += ((TVectorD *) (input->GetErrorOnBackground()->At(channel)))->operator[](bin) *
00304 toss[fgSystNames->BinarySearch((TObjString*) (((TObjArray *) (input->GetErrorNames()->At(channel)))->At(bin)))];
00305 }
00306 if ((serrf[channel] < -1.0) || (berrf[channel] < -0.9)) {
00307 retoss = kTRUE;
00308 continue;
00309 }
00310 }
00311 delete[]toss;
00312 } while (retoss);
00313
00314
00315 output->SetOwner();
00316 for (Int_t channel = 0; channel <= input->GetSignal()->GetLast();
00317 channel++) {
00318 TH1 *newsignal = (TH1*)(output->GetSignal()->At(channel));
00319 TH1 *oldsignal = (TH1*)(input->GetSignal()->At(channel));
00320 if(stat)
00321 for(int i=1; i<=newsignal->GetNbinsX(); i++)
00322 newsignal->SetBinContent(i,oldsignal->GetBinContent(i)+generator->Gaus(0,oldsignal->GetBinError(i)));
00323 else
00324 for(int i=1; i<=newsignal->GetNbinsX(); i++)
00325 newsignal->SetBinContent(i,oldsignal->GetBinContent(i));
00326 newsignal->Scale(1 + serrf[channel]);
00327 newsignal->SetDirectory(0);
00328 TH1 *newbackground = (TH1*)(output->GetBackground()->At(channel));
00329 TH1 *oldbackground = (TH1*)(input->GetBackground()->At(channel));
00330 if(stat)
00331 for(int i=1; i<=newbackground->GetNbinsX(); i++)
00332 newbackground->SetBinContent(i,oldbackground->GetBinContent(i)+generator->Gaus(0,oldbackground->GetBinError(i)));
00333 else
00334 for(int i=1; i<=newbackground->GetNbinsX(); i++)
00335 newbackground->SetBinContent(i,oldbackground->GetBinContent(i));
00336 newbackground->Scale(1 + berrf[channel]);
00337 newbackground->SetDirectory(0);
00338 }
00339 delete[] serrf;
00340 delete[] berrf;
00341 return 1;
00342 }
00343
00344 TConfidenceLevel *TLimit::ComputeLimit(TH1* s, TH1* b, TH1* d,
00345 Int_t nmc, bool stat,
00346 TRandom * generator)
00347 {
00348
00349
00350 TLimitDataSource* lds = new TLimitDataSource(s,b,d);
00351 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
00352 delete lds;
00353 return out;
00354 }
00355
00356 TConfidenceLevel *TLimit::ComputeLimit(TH1* s, TH1* b, TH1* d,
00357 TVectorD* se, TVectorD* be, TObjArray* l,
00358 Int_t nmc, bool stat,
00359 TRandom * generator)
00360 {
00361
00362
00363 TLimitDataSource* lds = new TLimitDataSource(s,b,d,se,be,l);
00364 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
00365 delete lds;
00366 return out;
00367 }
00368
00369 TConfidenceLevel *TLimit::ComputeLimit(Double_t s, Double_t b, Int_t d,
00370 Int_t nmc,
00371 bool stat,
00372 TRandom * generator)
00373 {
00374
00375
00376 TH1D* sh = new TH1D("__sh","__sh",1,0,2);
00377 sh->Fill(1,s);
00378 TH1D* bh = new TH1D("__bh","__bh",1,0,2);
00379 bh->Fill(1,b);
00380 TH1D* dh = new TH1D("__dh","__dh",1,0,2);
00381 dh->Fill(1,d);
00382 TLimitDataSource* lds = new TLimitDataSource(sh,bh,dh);
00383 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
00384 delete lds;
00385 delete sh;
00386 delete bh;
00387 delete dh;
00388 return out;
00389 }
00390
00391 TConfidenceLevel *TLimit::ComputeLimit(Double_t s, Double_t b, Int_t d,
00392 TVectorD* se, TVectorD* be, TObjArray* l,
00393 Int_t nmc,
00394 bool stat,
00395 TRandom * generator)
00396 {
00397
00398
00399 TH1D* sh = new TH1D("__sh","__sh",1,0,2);
00400 sh->Fill(1,s);
00401 TH1D* bh = new TH1D("__bh","__bh",1,0,2);
00402 bh->Fill(1,b);
00403 TH1D* dh = new TH1D("__dh","__dh",1,0,2);
00404 dh->Fill(1,d);
00405 TLimitDataSource* lds = new TLimitDataSource(sh,bh,dh,se,be,l);
00406 TConfidenceLevel* out = ComputeLimit(lds,nmc,stat,generator);
00407 delete lds;
00408 delete sh;
00409 delete bh;
00410 delete dh;
00411 return out;
00412 }
00413
00414 Double_t TLimit::LogLikelihood(Double_t s, Double_t b, Double_t b2, Double_t d)
00415 {
00416
00417
00418 return d*TMath::Log((s+b)/b2);
00419 }