rs401c_FeldmanCousins.C

Go to the documentation of this file.
00001 /////////////////////////////////////////////////////////////////////////
00002 //
00003 // Produces an interval on the mean signal in a number counting
00004 // experiment with known background using the Feldman-Cousins technique.
00005 // date Jan. 2009
00006 // updated June 2010
00007 //
00008 // Using the RooStats FeldmanCousins tool with 200 bins 
00009 // it takes 1 min and the interval is [0.2625, 10.6125] 
00010 // with a step size of 0.075.
00011 // The interval in Feldman & Cousins's original paper is [.29, 10.81]
00012 //  Phys.Rev.D57:3873-3889,1998. 
00013 /////////////////////////////////////////////////////////////////////////
00014 
00015 #include "RooGlobalFunc.h"
00016 #include "RooStats/ConfInterval.h"
00017 #include "RooStats/PointSetInterval.h"
00018 #include "RooStats/ConfidenceBelt.h"
00019 #include "RooStats/FeldmanCousins.h"
00020 #include "RooStats/ModelConfig.h"
00021 
00022 #include "RooWorkspace.h"
00023 #include "RooDataSet.h"
00024 #include "RooRealVar.h"
00025 #include "RooConstVar.h"
00026 #include "RooAddition.h"
00027 
00028 #include "RooDataHist.h"
00029 
00030 #include "RooPoisson.h"
00031 #include "RooPlot.h"
00032 
00033 #include "TCanvas.h"
00034 #include "TTree.h"
00035 #include "TH1F.h"
00036 #include "TMarker.h"
00037 #include "TStopwatch.h"
00038 
00039 #include <iostream>
00040 
00041 // use this order for safety on library loading
00042 using namespace RooFit ;
00043 using namespace RooStats ;
00044 
00045 
00046 void rs401c_FeldmanCousins()
00047 {
00048   // to time the macro... about 30 s
00049   TStopwatch t;
00050   t.Start();
00051 
00052   // make a simple model
00053   RooRealVar x("x","", 1,0,50);
00054   RooRealVar mu("mu","", 2.5,0, 15); // with a limit on mu>=0
00055   RooConstVar b("b","", 3.);
00056   RooAddition mean("mean","",RooArgList(mu,b));
00057   RooPoisson pois("pois", "", x, mean);
00058   RooArgSet parameters(mu);
00059 
00060   // create a toy dataset
00061   RooDataSet* data = pois.generate(RooArgSet(x), 1);
00062   data->Print("v");
00063   
00064   TCanvas* dataCanvas = new TCanvas("dataCanvas");
00065   RooPlot* frame = x.frame();
00066   data->plotOn(frame);
00067   frame->Draw();
00068   dataCanvas->Update();
00069 
00070   RooWorkspace* w = new RooWorkspace();
00071   ModelConfig modelConfig("poissonProblem",w);
00072   modelConfig.SetPdf(pois);
00073   modelConfig.SetParametersOfInterest(parameters);
00074   modelConfig.SetObservables(RooArgSet(x));
00075   w->Print();
00076 
00077   //////// show use of Feldman-Cousins
00078   RooStats::FeldmanCousins fc(*data,modelConfig);
00079   fc.SetTestSize(.05); // set size of test
00080   fc.UseAdaptiveSampling(true);
00081   fc.FluctuateNumDataEntries(false); // number counting analysis: dataset always has 1 entry with N events observed
00082   fc.SetNBins(100); // number of points to test per parameter
00083 
00084   // use the Feldman-Cousins tool
00085   PointSetInterval* interval = (PointSetInterval*)fc.GetInterval();
00086 
00087   // make a canvas for plots
00088   TCanvas* intervalCanvas =  new TCanvas("intervalCanvas");
00089   
00090   std::cout << "is this point in the interval? " << 
00091     interval->IsInInterval(parameters) << std::endl;
00092 
00093   std::cout << "interval is ["<<  
00094     interval->LowerLimit(mu)  << ", "  <<
00095     interval->UpperLimit(mu) << "]" << endl;
00096 
00097   // using 200 bins it takes 1 min and the answer is
00098   // interval is [0.2625, 10.6125] with a step size of .075
00099   // The interval in Feldman & Cousins's original paper is [.29, 10.81]
00100   //  Phys.Rev.D57:3873-3889,1998. 
00101 
00102   // No dedicated plotting class yet, so do it by hand:
00103 
00104   RooDataHist* parameterScan = (RooDataHist*) fc.GetPointsToScan();
00105   TH1F* hist = (TH1F*) parameterScan->createHistogram("mu",30);
00106   hist->Draw();
00107 
00108  
00109   RooArgSet* tmpPoint;
00110   // loop over points to test
00111   for(Int_t i=0; i<parameterScan->numEntries(); ++i){
00112     //    cout << "on parameter point " << i << " out of " << parameterScan->numEntries() << endl;
00113      // get a parameter point from the list of points to test.
00114     tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
00115 
00116     TMarker* mark = new TMarker(tmpPoint->getRealValue("mu"), 1, 25);
00117     if (interval->IsInInterval( *tmpPoint ) ) 
00118       mark->SetMarkerColor(kBlue);
00119     else
00120       mark->SetMarkerColor(kRed);
00121 
00122     mark->Draw("s");
00123     //delete tmpPoint;
00124     //    delete mark;
00125   }
00126   t.Stop();
00127   t.Print();
00128     
00129 
00130 }

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