void upper_limit_Bayesian_MCMC(Model* model,double confidence,int Niters){
  cout<<"///////////////////////////////////////////////////////////////////////////////////////////"<<endl;
  cout<<"Calculating upper limit with the MCMC method"<<endl;
  cout<<"///////////////////////////////////////////////////////////////////////////////////////////"<<endl;
  
  RooWorkspace* wspace = new RooWorkspace("wspace");
  ModelConfig* modelConfig = new ModelConfig("bayes");
  modelConfig->SetWorkspace(*wspace);
  modelConfig->SetPdf(*model->get_complete_likelihood());
  modelConfig->SetPriorPdf(*model->get_POI_prior());
  modelConfig->SetParametersOfInterest(*model->get_POI_set());
  //  modelConfig->SetNuisanceParameters();
  modelConfig->SetNuisanceParameters(*model->get_nuisance_set());


  //configure the calculator
  //model->Print();

  cout<<" POI size "<<model->get_POI_set()->getSize()<<endl; 
  RooRealVar* firstPOI = (RooRealVar*) modelConfig->GetParametersOfInterest()->first();


  MCMCCalculator mcmccalc(*model->get_data(), *modelConfig );
  mcmccalc.SetTestSize(1-confidence);
  mcmccalc.SetLeftSideTailFraction(0);
  mcmccalc.SetNumIters(Niters);
  MCMCInterval* interval = mcmccalc.GetInterval();
  double ul = interval->UpperLimit(*firstPOI);
  cout<<"UpperLimit: "<<ul<<endl;

}
Beispiel #2
0
double Tprime::printMcmcUpperLimit( double peak, std::string filename ) {
    //
    // print out the upper limit on the first Parameter of Interest
    //

    RooStats::ModelConfig * _mc = (RooStats::ModelConfig *)pWs->genobj("ModelConfig");

    RooRealVar * firstPOI = (RooRealVar*) _mc->GetParametersOfInterest()->first();
    double _limit = mcInt->UpperLimit(*firstPOI);
    cout << "\n95% upper limit on " <<firstPOI->GetName()<<" is : "<<
         _limit <<endl;

    if (filename.size()!=0) {

        std::ofstream aFile;

        // append to file if exists
        aFile.open(filename.c_str(), std::ios_base::app);

        char buf[1024];
        sprintf(buf, "%7.1f   %7.6f", peak, _limit);

        aFile << buf << std::endl;

        // close outfile here so it is safe even if subsequent iterations crash
        aFile.close();

    }

    return _limit;
}
Beispiel #3
0
double TwoBody::printMcmcUpperLimit( double peak, ModelConfig &_mc, std::string filename ){
  //
  // print out the upper limit on the first Parameter of Interest
  //

  std::string _legend = "[TwoBody::printMcmcUpperLimit]: ";

  char buf[1024];

  double _limit = numeric_limits<double>::max();

  if (mcInt){
    //mc.SetWorkspace(*ws);
    RooRealVar * firstPOI = (RooRealVar*) _mc.GetParametersOfInterest()->first();
    _limit = mcInt->UpperLimit(*firstPOI);
    std::cout << "\n95% upper limit on " << firstPOI->GetName() << " is : " <<
      _limit << endl;
    
    if (bMcmcConverged){
      sprintf(buf, "%7.1f   %7.6f", peak, _limit);
    }
    else{
      sprintf(buf, "# %7.1f   %7.6f  # MCMC did not converge", peak, _limit);
    }

  }
  else{
    sprintf(buf, "# MCMC did not converge");
  }

  if (filename.size()!=0){
    
    std::ofstream aFile;
    
    // append to file if exists
    aFile.open(filename.c_str(), std::ios_base::app);
    
    aFile << buf << std::endl;
    
    // close outfile here so it is safe even if subsequent iterations crash
    aFile.close();
    
  }
    
  return _limit;
}
Beispiel #4
0
MCMCInterval * TwoBody::GetMcmcInterval(ModelConfig _mc,
					double conf_level,
					int n_iter,
					int n_burn,
					double left_side_tail_fraction,
					int n_bins){
  // use MCMCCalculator  (takes about 1 min)
  // Want an efficient proposal function, so derive it from covariance
  // matrix of fit

  std::string legend = "[TwoBody::GetMcmcInterval]: ";

  // FIXME: testing: this proposal function seems fairly robust
  SequentialProposal sp(0.5);

  MCMCCalculator mcmc( *data, _mc );
  mcmc.SetConfidenceLevel(conf_level);
  mcmc.SetNumIters(n_iter);          // Metropolis-Hastings algorithm iterations

  mcmc.SetProposalFunction(sp);

  mcmc.SetNumBurnInSteps(n_burn); // first N steps to be ignored as burn-in
  mcmc.SetLeftSideTailFraction(left_side_tail_fraction);
  mcmc.SetNumBins(n_bins);
  
  // FIXME: testing good initial values - don't seem to do anything different
  //ws->var("ratio")->setVal(0.01);
  //ws->var("beta_nsig_dielectron")->setRange(-3.0, 3.0);
  //ws->var("beta_nbkg_dielectron")->setRange(-3.0, 3.0);
  //ws->var("beta_mass_dielectron")->setRange(-3.0, 3.0);

//mcInt = mcmc.GetInterval();
  try {
    delete mcInt;
    mcInt = mcmc.GetInterval();
  } catch ( std::length_error &ex) {
    mcInt = 0;
  }
  cerr<<legend<<endl;

  // check if limit makes sense
  bMcmcConverged = false; // default
  if (mcInt){
    RooRealVar * p_first_poi = (RooRealVar*) _mc.GetParametersOfInterest()->first();
    double poi_limit = mcInt->UpperLimit(*p_first_poi);
    double u_poi_min  = p_first_poi->getMin();
    double u_poi_max  = p_first_poi->getMax();
    double u_poi_gap = (u_poi_max-poi_limit)/(u_poi_max-u_poi_min);
    std::cout << legend << "POI upper limit: " << poi_limit << std::endl;
    std::cout << legend << "POI range: [" << u_poi_min << ", " << u_poi_max << "]" << std::endl;
    std::cout << legend << "POI upper gap (fraction of range): " << u_poi_gap << std::endl;
    if (u_poi_gap<0.2){
      std::cout << legend 
		<< "POI limit too close to the upper boundary, MCMC probably failed!!!" << std::endl;
      std::cout << legend
		<< "returning interval and setting fail flag" << std::endl;
      bMcmcConverged = false;
    }
    else{
      bMcmcConverged = true;
    }
  }
  else std::cout << "No interval found!" << std::endl;
  
  return mcInt;
}
Beispiel #5
0
//
// calculation of the limit: assumes that wspace is set up and observations
//   contained in data
//
MyLimit computeLimit (RooWorkspace* wspace, RooDataSet* data, StatMethod method, bool draw) {

  // let's time this challenging example
  TStopwatch t;

  //
  // get nominal signal
  //
  RooRealVar exp_sig(*wspace->var("s"));
  double exp_sig_val = exp_sig.getVal();
  std::cout << "exp_sig = " << exp_sig_val << std::endl;
  
  /////////////////////////////////////////////////////
  // Now the statistical tests
  // model config
  std::cout << wspace->pdf("model") << " "
	    << wspace->pdf("prior") << " "
	    << wspace->set("poi") << " "
	    << wspace->set("nuis") << std::endl;
  ModelConfig modelConfig("RA4abcd");
  modelConfig.SetWorkspace(*wspace);
  modelConfig.SetPdf(*wspace->pdf("model"));
  modelConfig.SetPriorPdf(*wspace->pdf("prior"));
  modelConfig.SetParametersOfInterest(*wspace->set("poi"));
  modelConfig.SetNuisanceParameters(*wspace->set("nuis"));


  //////////////////////////////////////////////////
  // If you want to see the covariance matrix uncomment
  // wspace->pdf("model")->fitTo(*data);

  // use ProfileLikelihood
  if ( method == ProfileLikelihoodMethod ) {
    ProfileLikelihoodCalculator plc(*data, modelConfig);
    plc.SetConfidenceLevel(0.95);
    RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
    RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
    LikelihoodInterval* plInt = plc.GetInterval();
    RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
    plInt->LowerLimit( *wspace->var("s") ); // get ugly print out of the way. Fix.
    // RooMsgService::instance().setGlobalKillBelow(RooFit::DEBUG);
    if ( draw ) {
      TCanvas* c = new TCanvas("ProfileLikelihood");
      LikelihoodIntervalPlot* lrplot = new LikelihoodIntervalPlot(plInt);
      lrplot->Draw();
    }
//     RooMsgService::instance().setGlobalKillBelow(msglevel);
    double lowLim = plInt->LowerLimit(*wspace->var("s"));
    double uppLim = plInt->UpperLimit(*wspace->var("s"));
//     double exp_sig_val = wspace->var("s")->getVal();
//     double exp_sig_val = exp_sig.getVal();
    cout << "Profile Likelihood interval on s = [" << 
      lowLim << ", " <<
      uppLim << "]" << " " << exp_sig_val << endl; 
//     MyLimit result(plInt->IsInInterval(exp_sig),
    MyLimit result(exp_sig_val>lowLim&&exp_sig_val<uppLim,lowLim,uppLim);
    // std::cout << "isIn " << result << std::endl;
    delete plInt;
//     delete modelConfig;
    return result;
  }

  // use FeldmaCousins (takes ~20 min)  
  if ( method == FeldmanCousinsMethod ) {
    FeldmanCousins fc(*data, modelConfig);
    fc.SetConfidenceLevel(0.95);
    //number counting: dataset always has 1 entry with N events observed
    fc.FluctuateNumDataEntries(false); 
    fc.UseAdaptiveSampling(true);
    fc.SetNBins(100);
    PointSetInterval* fcInt = NULL;
    fcInt = (PointSetInterval*) fc.GetInterval(); // fix cast
    double lowLim = fcInt->LowerLimit(*wspace->var("s"));
    double uppLim = fcInt->UpperLimit(*wspace->var("s"));
//     double exp_sig_val = wspace->var("s")->getVal();
    cout << "Feldman Cousins interval on s = [" << lowLim << " " << uppLim << endl;
    // std::cout << "isIn " << result << std::endl;
    MyLimit result(exp_sig_val>lowLim&&exp_sig_val<uppLim,
		   fcInt->LowerLimit(*wspace->var("s")),fcInt->UpperLimit(*wspace->var("s")));
    delete fcInt;
    return result;
  }


  // use BayesianCalculator (only 1-d parameter of interest, slow for this problem)  
  if ( method == BayesianMethod ) {
    BayesianCalculator bc(*data, modelConfig);
    bc.SetConfidenceLevel(0.95);
    bc.SetLeftSideTailFraction(0.5);
    SimpleInterval* bInt = NULL;
    if( wspace->set("poi")->getSize() == 1)   {
      bInt = bc.GetInterval();
      if ( draw ) {
	TCanvas* c = new TCanvas("Bayesian");
	// the plot takes a long time and print lots of error
	// using a scan it is better
	bc.SetScanOfPosterior(50);
	RooPlot* bplot = bc.GetPosteriorPlot();
	bplot->Draw();
      }
      cout << "Bayesian interval on s = [" << 
	bInt->LowerLimit( ) << ", " <<
	bInt->UpperLimit( ) << "]" << endl;
      // std::cout << "isIn " << result << std::endl;
      MyLimit result(bInt->IsInInterval(exp_sig),
		     bInt->LowerLimit(),bInt->UpperLimit());
      delete bInt;
      return result;
    } else {
    cout << "Bayesian Calc. only supports on parameter of interest" << endl;
    return MyLimit();
    }
  }


  // use MCMCCalculator  (takes about 1 min)
  // Want an efficient proposal function, so derive it from covariance
  // matrix of fit
  if ( method == MCMCMethod ) {
    RooFitResult* fit = wspace->pdf("model")->fitTo(*data,Save());
    ProposalHelper ph;
    ph.SetVariables((RooArgSet&)fit->floatParsFinal());
    ph.SetCovMatrix(fit->covarianceMatrix());
    ph.SetUpdateProposalParameters(kTRUE); // auto-create mean vars and add mappings
    ph.SetCacheSize(100);
    ProposalFunction* pf = ph.GetProposalFunction();
    
    MCMCCalculator mc(*data, modelConfig);
    mc.SetConfidenceLevel(0.95);
    mc.SetProposalFunction(*pf);
    mc.SetNumBurnInSteps(100); // first N steps to be ignored as burn-in
    mc.SetNumIters(100000);
    mc.SetLeftSideTailFraction(0.5); // make a central interval
    MCMCInterval* mcInt = NULL;
    mcInt = mc.GetInterval();
    MCMCIntervalPlot mcPlot(*mcInt); 
    mcPlot.Draw();
    cout << "MCMC interval on s = [" << 
      mcInt->LowerLimit(*wspace->var("s") ) << ", " <<
      mcInt->UpperLimit(*wspace->var("s") ) << "]" << endl;
    // std::cout << "isIn " << result << std::endl;
    MyLimit result(mcInt->IsInInterval(exp_sig),
		   mcInt->LowerLimit(*wspace->var("s")),mcInt->UpperLimit(*wspace->var("s")));
    delete mcInt;
    return result;
  }
  

  t.Print();

//   delete modelConfig;
  return MyLimit();

}
Beispiel #6
0
void IntervalExamples()
{

   // Time this macro
   TStopwatch t;
   t.Start();


   // set RooFit random seed for reproducible results
   RooRandom::randomGenerator()->SetSeed(3001);

   // make a simple model via the workspace factory
   RooWorkspace* wspace = new RooWorkspace();
   wspace->factory("Gaussian::normal(x[-10,10],mu[-1,1],sigma[1])");
   wspace->defineSet("poi","mu");
   wspace->defineSet("obs","x");

   // specify components of model for statistical tools
   ModelConfig* modelConfig = new ModelConfig("Example G(x|mu,1)");
   modelConfig->SetWorkspace(*wspace);
   modelConfig->SetPdf( *wspace->pdf("normal") );
   modelConfig->SetParametersOfInterest( *wspace->set("poi") );
   modelConfig->SetObservables( *wspace->set("obs") );

   // create a toy dataset
   RooDataSet* data = wspace->pdf("normal")->generate(*wspace->set("obs"),100);
   data->Print();

   // for convenience later on
   RooRealVar* x = wspace->var("x");
   RooRealVar* mu = wspace->var("mu");

   // set confidence level
   double confidenceLevel = 0.95;

   // example use profile likelihood calculator
   ProfileLikelihoodCalculator plc(*data, *modelConfig);
   plc.SetConfidenceLevel( confidenceLevel);
   LikelihoodInterval* plInt = plc.GetInterval();

   // example use of Feldman-Cousins
   FeldmanCousins fc(*data, *modelConfig);
   fc.SetConfidenceLevel( confidenceLevel);
   fc.SetNBins(100); // number of points to test per parameter
   fc.UseAdaptiveSampling(true); // make it go faster

   // Here, we consider only ensembles with 100 events
   // The PDF could be extended and this could be removed
   fc.FluctuateNumDataEntries(false);

   // Proof
   //  ProofConfig pc(*wspace, 4, "workers=4", kFALSE);    // proof-lite
   //ProofConfig pc(w, 8, "localhost");    // proof cluster at "localhost"
   //  ToyMCSampler* toymcsampler = (ToyMCSampler*) fc.GetTestStatSampler();
   //  toymcsampler->SetProofConfig(&pc);     // enable proof

   PointSetInterval* interval = (PointSetInterval*) fc.GetInterval();


   // example use of BayesianCalculator
   // now we also need to specify a prior in the ModelConfig
   wspace->factory("Uniform::prior(mu)");
   modelConfig->SetPriorPdf(*wspace->pdf("prior"));

   // example usage of BayesianCalculator
   BayesianCalculator bc(*data, *modelConfig);
   bc.SetConfidenceLevel( confidenceLevel);
   SimpleInterval* bcInt = bc.GetInterval();

   // example use of MCMCInterval
   MCMCCalculator mc(*data, *modelConfig);
   mc.SetConfidenceLevel( confidenceLevel);
   // special options
   mc.SetNumBins(200);        // bins used internally for representing posterior
   mc.SetNumBurnInSteps(500); // first N steps to be ignored as burn-in
   mc.SetNumIters(100000);    // how long to run chain
   mc.SetLeftSideTailFraction(0.5); // for central interval
   MCMCInterval* mcInt = mc.GetInterval();

   // for this example we know the expected intervals
   double expectedLL = data->mean(*x)
      + ROOT::Math::normal_quantile(  (1-confidenceLevel)/2,1)
      / sqrt(data->numEntries());
   double expectedUL = data->mean(*x)
      + ROOT::Math::normal_quantile_c((1-confidenceLevel)/2,1)
      / sqrt(data->numEntries()) ;

   // Use the intervals
   std::cout << "expected interval is [" <<
      expectedLL << ", " <<
      expectedUL << "]" << endl;

   cout << "plc interval is [" <<
      plInt->LowerLimit(*mu) << ", " <<
      plInt->UpperLimit(*mu) << "]" << endl;

   std::cout << "fc interval is ["<<
      interval->LowerLimit(*mu) << " , "  <<
      interval->UpperLimit(*mu) << "]" << endl;

   cout << "bc interval is [" <<
      bcInt->LowerLimit() << ", " <<
      bcInt->UpperLimit() << "]" << endl;

   cout << "mc interval is [" <<
      mcInt->LowerLimit(*mu) << ", " <<
      mcInt->UpperLimit(*mu) << "]" << endl;

   mu->setVal(0);
   cout << "is mu=0 in the interval? " <<
      plInt->IsInInterval(RooArgSet(*mu)) << endl;


   // make a reasonable style
   gStyle->SetCanvasColor(0);
   gStyle->SetCanvasBorderMode(0);
   gStyle->SetPadBorderMode(0);
   gStyle->SetPadColor(0);
   gStyle->SetCanvasColor(0);
   gStyle->SetTitleFillColor(0);
   gStyle->SetFillColor(0);
   gStyle->SetFrameFillColor(0);
   gStyle->SetStatColor(0);


   // some plots
   TCanvas* canvas = new TCanvas("canvas");
   canvas->Divide(2,2);

   // plot the data
   canvas->cd(1);
   RooPlot* frame = x->frame();
   data->plotOn(frame);
   data->statOn(frame);
   frame->Draw();

   // plot the profile likelihood
   canvas->cd(2);
   LikelihoodIntervalPlot plot(plInt);
   plot.Draw();

   // plot the MCMC interval
   canvas->cd(3);
   MCMCIntervalPlot* mcPlot = new MCMCIntervalPlot(*mcInt);
   mcPlot->SetLineColor(kGreen);
   mcPlot->SetLineWidth(2);
   mcPlot->Draw();

   canvas->cd(4);
   RooPlot * bcPlot = bc.GetPosteriorPlot();
   bcPlot->Draw();

   canvas->Update();

   t.Stop();
   t.Print();

}
Beispiel #7
0
Int_t Tprime::RunMcmc( std::string channel, // ejets, mujets, combined
                       std::string mode,    // observed, expected
                       double peak,         // resonance mass
                       std::string suffix,  // suffix for output file names
                       Int_t ntoys,         // number of pseudoexperiments for expected limit
                       Int_t mcmc_iter,     // number of MCMC iterations
                       Int_t mcmc_burnin,   // number of MCMC burn in steps to be discarded
                       std::string inputdir // input dir name
                     ) {
    //
    // Bayesian MCMC calculation
    //

    std::string legend = "[tprime::RunMcmc()]: ";

    // print out inputs
    std::cout << legend << std::endl;
    std::cout << legend << "Input parameters specified. Some of them are not used and defaults are entered" << std::endl;
    std::cout << legend << "------------------------------" << std::endl;
    std::cout << legend << "channel: " << channel << std::endl;
    std::cout << legend << "mode: " << mode << std::endl;
    std::cout << legend << "input directory: " << inputdir << std::endl;
    std::cout << legend << "resonance peak mass: " << peak << std::endl;
    std::cout << legend << "suffix: ->" << suffix << "<-" << std::endl;
    std::cout << legend << "number of pseudo-experiments: "<< ntoys << std::endl;
    std::cout << legend << std::endl;
    std::cout << legend << "Bayesian MCMC parameters" << std::endl;
    std::cout << legend << "------------------------------" << std::endl;
    std::cout << legend << "number of iterations: " << mcmc_iter << std::endl;
    std::cout << legend << "number of burn-in steps to discard: " << mcmc_burnin << std::endl;
    std::cout << legend << std::endl;

    // compose the workspace file name
    char buf[1024];
    sprintf(buf, "%sresults_%04.0f/tprime_%s_tprimeCrossSection_model.root", inputdir.c_str(), peak, channel.c_str());
    std::string _file = buf;
    std::cout << legend << "guessed name of the file with the workspace: >" << _file << "<" << std::endl;

    //load workspace
    LoadWorkspace(_file, channel);

    // change POI range
    double poiUpper = GetPoiUpper(channel, peak);
    std::cout << legend << "setting POI range to [0; " << poiUpper << "]" << std::endl;
    pWs->var("xsec")->setRange(0.0, poiUpper);

    // timer
    TStopwatch t;
    t.Start();

    int pe_counter = 0;
    std::vector<Double_t> _limits;
    while (pe_counter < ntoys) {

        // for mass limit, add k-factor systematics to the nsig systematics
        // FIXME: this is a correlated part of the uncertainty!!!
        //  - different uncertainties for graviton and Z' models
        if ( mode.find("mass_") != std::string::npos ) {

            std::cout << legend << std::endl;
            std::cout << legend << "this a mass limit calculation," << std::endl;
            std::cout << legend << "I would add k-factor uncertainty to the nsig uncertainty" << std::endl;
            std::cout << legend << "I would do it " << ntoys << " times, so one can average. " << pe_counter+1 << " of " << ntoys << std::endl;
            std::cout << legend << "Not implemented yet " << std::endl;
            std::cout << legend << std::endl;

            //Double_t kfactor_err = GetKfactorUncertainty(peak, mode);

            //double nsig_kappa = ws->var("nsig_kappa_dimuon")->getVal();
            //nsig_kappa = sqrt(nsig_kappa*nsig_kappa + kfactor_err*kfactor_err);
            //ws->var("nsig_kappa_dimuon")->setVal(nsig_kappa);

            //ntoys = 1;

        }

        else if ( mode.find("expected") != std::string::npos ) {

            std::cout << legend << std::endl;
            std::cout << legend << "this is pseudoexperiment " << pe_counter+1 << " of " << ntoys << std::endl;
            std::cout << legend << "for the expected limit estimate" << std::endl;
            std::cout << legend << std::endl;

            // prepare PE data
            GetPseudoData();

        }

        else { //  "regular" observed limit

            std::cout << legend << std::endl;
            std::cout << legend << "calculating an observed limit..." << std::endl;
            std::cout << legend << "I will do it " << ntoys << " times, so one can average. " << pe_counter+1 << " of " << ntoys << std::endl;
            std::cout << legend << std::endl;

            GetWorkspaceData("obsData");

            //ntoys = 1;
        }

        mcInt = GetMcmcInterval(0.95,        // conf level
                                mcmc_iter,   // number of iterations
                                mcmc_burnin, // number of burn-in to discard
                                0.0,         // left side tail fraction, 0 for upper limit
                                100);        // number of bins in posterior, only for plotting

        ++pe_counter;

        if (!mcInt) {
            continue;
        }
        else {

            std::string _outfile = "tprime_"+channel+"_xsec_mcmc_limit_" + suffix + ".ascii";
            printMcmcUpperLimit( peak, _outfile );

            // limits for averaging/medianing
            RooStats::ModelConfig * pSbModel = GetModelConfig("ModelConfig");
            RooRealVar * firstPOI = (RooRealVar*) pSbModel->GetParametersOfInterest()->first();
            double _limit = mcInt->UpperLimit(*firstPOI);
            _limits.push_back(_limit);

        } // end of valid mcInt block

    } // end of while

    // write median limit to a file
    if (_limits.size() > 0) {
        Double_t _median_limit = TMath::Median(_limits.size(), &_limits[0]);
        std::vector<Double_t> _mass_limit;
        _mass_limit.push_back(peak);
        _mass_limit.push_back(_median_limit);
        std::string _outfile = "tprime_"+channel+"_xsec_mcmc_median_limit_" + suffix + ".ascii";
        PrintToFile(_outfile, _mass_limit, "#  mass         median_limit");
    }

    std::string _outfile = "tprime_"+channel+"_xsec_mcmc_posterior_" + suffix + ".pdf";
    makeMcmcPosteriorPlot( _outfile );

    // timer
    t.Print();

    return 0;
}
Beispiel #8
0
void rs101_limitexample()
{
    // --------------------------------------
    // An example of setting a limit in a number counting experiment with uncertainty on background and signal

    // to time the macro
    TStopwatch t;
    t.Start();

    // --------------------------------------
    // The Model building stage
    // --------------------------------------
    RooWorkspace* wspace = new RooWorkspace();
    wspace->factory("Poisson::countingModel(obs[150,0,300], sum(s[50,0,120]*ratioSigEff[1.,0,3.],b[100]*ratioBkgEff[1.,0.,3.]))"); // counting model
    //  wspace->factory("Gaussian::sigConstraint(ratioSigEff,1,0.05)"); // 5% signal efficiency uncertainty
    //  wspace->factory("Gaussian::bkgConstraint(ratioBkgEff,1,0.1)"); // 10% background efficiency uncertainty
    wspace->factory("Gaussian::sigConstraint(gSigEff[1,0,3],ratioSigEff,0.05)"); // 5% signal efficiency uncertainty
    wspace->factory("Gaussian::bkgConstraint(gSigBkg[1,0,3],ratioBkgEff,0.2)"); // 10% background efficiency uncertainty
    wspace->factory("PROD::modelWithConstraints(countingModel,sigConstraint,bkgConstraint)"); // product of terms
    wspace->Print();

    RooAbsPdf* modelWithConstraints = wspace->pdf("modelWithConstraints"); // get the model
    RooRealVar* obs = wspace->var("obs"); // get the observable
    RooRealVar* s = wspace->var("s"); // get the signal we care about
    RooRealVar* b = wspace->var("b"); // get the background and set it to a constant.  Uncertainty included in ratioBkgEff
    b->setConstant();

    RooRealVar* ratioSigEff = wspace->var("ratioSigEff"); // get uncertain parameter to constrain
    RooRealVar* ratioBkgEff = wspace->var("ratioBkgEff"); // get uncertain parameter to constrain
    RooArgSet constrainedParams(*ratioSigEff, *ratioBkgEff); // need to constrain these in the fit (should change default behavior)

    RooRealVar * gSigEff = wspace->var("gSigEff");     // global observables for signal efficiency
    RooRealVar * gSigBkg = wspace->var("gSigBkg");  // global obs for background efficiency
    gSigEff->setConstant();
    gSigBkg->setConstant();

    // Create an example dataset with 160 observed events
    obs->setVal(160.);
    RooDataSet* data = new RooDataSet("exampleData", "exampleData", RooArgSet(*obs));
    data->add(*obs);

    RooArgSet all(*s, *ratioBkgEff, *ratioSigEff);

    // not necessary
    modelWithConstraints->fitTo(*data, RooFit::Constrain(RooArgSet(*ratioSigEff, *ratioBkgEff)));

    // Now let's make some confidence intervals for s, our parameter of interest
    RooArgSet paramOfInterest(*s);

    ModelConfig modelConfig(wspace);
    modelConfig.SetPdf(*modelWithConstraints);
    modelConfig.SetParametersOfInterest(paramOfInterest);
    modelConfig.SetNuisanceParameters(constrainedParams);
    modelConfig.SetObservables(*obs);
    modelConfig.SetGlobalObservables( RooArgSet(*gSigEff,*gSigBkg));
    modelConfig.SetName("ModelConfig");
    wspace->import(modelConfig);
    wspace->import(*data);
    wspace->SetName("w");
    wspace->writeToFile("rs101_ws.root");



    // First, let's use a Calculator based on the Profile Likelihood Ratio
    //ProfileLikelihoodCalculator plc(*data, *modelWithConstraints, paramOfInterest);
    ProfileLikelihoodCalculator plc(*data, modelConfig);
    plc.SetTestSize(.05);
    ConfInterval* lrinterval = plc.GetInterval();  // that was easy.

    // Let's make a plot
    TCanvas* dataCanvas = new TCanvas("dataCanvas");
    dataCanvas->Divide(2,1);

    dataCanvas->cd(1);
    LikelihoodIntervalPlot plotInt((LikelihoodInterval*)lrinterval);
    plotInt.SetTitle("Profile Likelihood Ratio and Posterior for S");
    plotInt.Draw();

    // Second, use a Calculator based on the Feldman Cousins technique
    FeldmanCousins fc(*data, modelConfig);
    fc.UseAdaptiveSampling(true);
    fc.FluctuateNumDataEntries(false); // number counting analysis: dataset always has 1 entry with N events observed
    fc.SetNBins(100); // number of points to test per parameter
    fc.SetTestSize(.05);
    //  fc.SaveBeltToFile(true); // optional
    ConfInterval* fcint = NULL;
    fcint = fc.GetInterval();  // that was easy.

    RooFitResult* fit = modelWithConstraints->fitTo(*data, Save(true));

    // Third, use a Calculator based on Markov Chain monte carlo
    // Before configuring the calculator, let's make a ProposalFunction
    // that will achieve a high acceptance rate
    ProposalHelper ph;
    ph.SetVariables((RooArgSet&)fit->floatParsFinal());
    ph.SetCovMatrix(fit->covarianceMatrix());
    ph.SetUpdateProposalParameters(true);
    ph.SetCacheSize(100);
    ProposalFunction* pdfProp = ph.GetProposalFunction();  // that was easy

    MCMCCalculator mc(*data, modelConfig);
    mc.SetNumIters(20000); // steps to propose in the chain
    mc.SetTestSize(.05); // 95% CL
    mc.SetNumBurnInSteps(40); // ignore first N steps in chain as "burn in"
    mc.SetProposalFunction(*pdfProp);
    mc.SetLeftSideTailFraction(0.5);  // find a "central" interval
    MCMCInterval* mcInt = (MCMCInterval*)mc.GetInterval();  // that was easy


    // Get Lower and Upper limits from Profile Calculator
    cout << "Profile lower limit on s = " << ((LikelihoodInterval*) lrinterval)->LowerLimit(*s) << endl;
    cout << "Profile upper limit on s = " << ((LikelihoodInterval*) lrinterval)->UpperLimit(*s) << endl;

    // Get Lower and Upper limits from FeldmanCousins with profile construction
    if (fcint != NULL) {
        double fcul = ((PointSetInterval*) fcint)->UpperLimit(*s);
        double fcll = ((PointSetInterval*) fcint)->LowerLimit(*s);
        cout << "FC lower limit on s = " << fcll << endl;
        cout << "FC upper limit on s = " << fcul << endl;
        TLine* fcllLine = new TLine(fcll, 0, fcll, 1);
        TLine* fculLine = new TLine(fcul, 0, fcul, 1);
        fcllLine->SetLineColor(kRed);
        fculLine->SetLineColor(kRed);
        fcllLine->Draw("same");
        fculLine->Draw("same");
        dataCanvas->Update();
    }

    // Plot MCMC interval and print some statistics
    MCMCIntervalPlot mcPlot(*mcInt);
    mcPlot.SetLineColor(kMagenta);
    mcPlot.SetLineWidth(2);
    mcPlot.Draw("same");

    double mcul = mcInt->UpperLimit(*s);
    double mcll = mcInt->LowerLimit(*s);
    cout << "MCMC lower limit on s = " << mcll << endl;
    cout << "MCMC upper limit on s = " << mcul << endl;
    cout << "MCMC Actual confidence level: "
         << mcInt->GetActualConfidenceLevel() << endl;

    // 3-d plot of the parameter points
    dataCanvas->cd(2);
    // also plot the points in the markov chain
    RooDataSet * chainData = mcInt->GetChainAsDataSet();

    assert(chainData);
    std::cout << "plotting the chain data - nentries = " << chainData->numEntries() << std::endl;
    TTree* chain =  RooStats::GetAsTTree("chainTreeData","chainTreeData",*chainData);
    assert(chain);
    chain->SetMarkerStyle(6);
    chain->SetMarkerColor(kRed);

    chain->Draw("s:ratioSigEff:ratioBkgEff","nll_MarkovChain_local_","box"); // 3-d box proportional to posterior

    // the points used in the profile construction
    RooDataSet * parScanData = (RooDataSet*) fc.GetPointsToScan();
    assert(parScanData);
    std::cout << "plotting the scanned points used in the frequentist construction - npoints = " << parScanData->numEntries() << std::endl;
    // getting the tree and drawing it -crashes (very strange....);
    // TTree* parameterScan =  RooStats::GetAsTTree("parScanTreeData","parScanTreeData",*parScanData);
    // assert(parameterScan);
    // parameterScan->Draw("s:ratioSigEff:ratioBkgEff","","goff");
    TGraph2D *gr = new TGraph2D(parScanData->numEntries());
    for (int ievt = 0; ievt < parScanData->numEntries(); ++ievt) {
        const RooArgSet * evt = parScanData->get(ievt);
        double x = evt->getRealValue("ratioBkgEff");
        double y = evt->getRealValue("ratioSigEff");
        double z = evt->getRealValue("s");
        gr->SetPoint(ievt, x,y,z);
        // std::cout << ievt << "  " << x << "  " << y << "  " << z << std::endl;
    }
    gr->SetMarkerStyle(24);
    gr->Draw("P SAME");


    delete wspace;
    delete lrinterval;
    delete mcInt;
    delete fcint;
    delete data;

    // print timing info
    t.Stop();
    t.Print();
}