limit.C

Go to the documentation of this file.
00001 #include <iostream>
00002 #include "TH1.h"
00003 #include "THStack.h"
00004 #include "TCanvas.h"
00005 #include "TFrame.h"
00006 #include "TRandom2.h"
00007 #include "TSystem.h"
00008 #include "TVector.h"
00009 #include "TObjArray.h"
00010 #include "TLimit.h"
00011 #include "TLimitDataSource.h"
00012 #include "TConfidenceLevel.h"
00013 
00014 using std::cout; 
00015 using std::endl;
00016 
00017 void limit() {
00018 //This program demonstrates the computation of 95 % C.L. limits.
00019 //It uses a set of randomly created histograms.
00020 //
00021 //Author:  Christophe.Delaere@cern.ch on 21.08.02
00022 
00023 // Create a new canvas.
00024   TCanvas *c1 = new TCanvas("c1","Dynamic Filling Example",200,10,700,500);
00025   c1->SetFillColor(42);
00026           
00027 // Create some histograms
00028   TH1D* background = new TH1D("background","The expected background",30,-4,4);
00029   TH1D* signal     = new TH1D("signal","the expected signal",30,-4,4);
00030   TH1D* data       = new TH1D("data","some fake data points",30,-4,4);
00031   background->SetFillColor(48);
00032   signal->SetFillColor(41);
00033   data->SetMarkerStyle(21);
00034   data->SetMarkerColor(kBlue);
00035   background->Sumw2(); // needed for stat uncertainty
00036   signal->Sumw2(); // needed for stat uncertainty
00037   
00038 // Fill histograms randomly
00039   TRandom2 r;
00040   Float_t bg,sig,dt;
00041   for (Int_t i = 0; i < 25000; i++) {
00042      bg  = r.Gaus(0,1);
00043      sig = r.Gaus(1,.2);
00044      background->Fill(bg,0.02);
00045      signal->Fill(sig,0.001);
00046   }
00047   for (Int_t i = 0; i < 500; i++) {
00048      dt = r.Gaus(0,1);
00049      data->Fill(dt);
00050   }
00051   THStack *hs = new THStack("hs","Signal and background compared to data...");
00052   hs->Add(background);
00053   hs->Add(signal);
00054   hs->Draw("hist");
00055   data->Draw("PE1,Same");
00056   c1->Modified();
00057   c1->Update();
00058   c1->GetFrame()->SetFillColor(21);
00059   c1->GetFrame()->SetBorderSize(6);
00060   c1->GetFrame()->SetBorderMode(-1);
00061   c1->Modified();
00062   c1->Update();
00063   gSystem->ProcessEvents();
00064 
00065 // Compute the limits
00066   cout << "Computing limits... " << endl;
00067   TLimitDataSource* mydatasource = new TLimitDataSource(signal,background,data);
00068   TConfidenceLevel *myconfidence = TLimit::ComputeLimit(mydatasource,50000);
00069   cout << "CLs    : "   << myconfidence->CLs()  << endl;
00070   cout << "CLsb   : "   << myconfidence->CLsb() << endl;
00071   cout << "CLb    : "   << myconfidence->CLb()  << endl;
00072   cout << "< CLs >  : " << myconfidence->GetExpectedCLs_b()  << endl;
00073   cout << "< CLsb > : " << myconfidence->GetExpectedCLsb_b() << endl;
00074   cout << "< CLb >  : " << myconfidence->GetExpectedCLb_b()  << endl;
00075 
00076 // Add stat uncertainty
00077   cout << endl << "Computing limits with stat systematics... " << endl;
00078   TConfidenceLevel *mystatconfidence = TLimit::ComputeLimit(mydatasource,50000,true);
00079   cout << "CLs    : "   << mystatconfidence->CLs()  << endl;
00080   cout << "CLsb   : "   << mystatconfidence->CLsb() << endl;
00081   cout << "CLb    : "   << mystatconfidence->CLb()  << endl;
00082   cout << "< CLs >  : " << mystatconfidence->GetExpectedCLs_b()  << endl;
00083   cout << "< CLsb > : " << mystatconfidence->GetExpectedCLsb_b() << endl;
00084   cout << "< CLb >  : " << mystatconfidence->GetExpectedCLb_b()  << endl;
00085 
00086 // Add some systematics
00087   cout << endl << "Computing limits with systematics... " << endl;
00088   TVectorD errorb(2);
00089   TVectorD errors(2);
00090   TObjArray* names = new TObjArray();
00091   TObjString name1("bg uncertainty");
00092   TObjString name2("sig uncertainty");
00093   names->AddLast(&name1);
00094   names->AddLast(&name2);
00095   errorb[0]=0.05; // error source 1: 5%
00096   errorb[1]=0;    // error source 2: 0%
00097   errors[0]=0;    // error source 1: 0%
00098   errors[1]=0.01; // error source 2: 1%
00099   TLimitDataSource* mynewdatasource  = new TLimitDataSource();
00100   mynewdatasource->AddChannel(signal,background,data,&errors,&errorb,names);
00101   TConfidenceLevel *mynewconfidence = TLimit::ComputeLimit(mynewdatasource,50000,true);
00102   cout << "CLs    : " << mynewconfidence->CLs()  << endl;
00103   cout << "CLsb   : " << mynewconfidence->CLsb() << endl;
00104   cout << "CLb    : " << mynewconfidence->CLb()  << endl;
00105   cout << "< CLs >  : " << mynewconfidence->GetExpectedCLs_b()  << endl;
00106   cout << "< CLsb > : " << mynewconfidence->GetExpectedCLsb_b() << endl;
00107   cout << "< CLb >  : " << mynewconfidence->GetExpectedCLb_b()  << endl;
00108 
00109 // show canonical -2lnQ plots in a new canvas
00110 // - The histogram of -2lnQ for background hypothesis (full)
00111 // - The histogram of -2lnQ for signal and background hypothesis (dashed)
00112   TCanvas *c2 = new TCanvas("c2");
00113   myconfidence->Draw();
00114   
00115 // clean up (except histograms and canvas)
00116   delete myconfidence;
00117   delete mydatasource;
00118   delete mystatconfidence;
00119   delete mynewconfidence;
00120   delete mynewdatasource;
00121 }
00122   

Generated on Tue Jul 5 15:44:49 2011 for ROOT_528-00b_version by  doxygen 1.5.1