Extended likelihood PDFs
// extended likelihood use
intro7()
{
// Build regular Gaussian PDF
RooRealVar x("x","x",-10,10) ;
RooRealVar mean("mean","mean of gaussian",-3,-10,10) ;
RooRealVar sigma("sigma","width of gaussian",1,0.1,5) ;
RooGaussian gauss("gauss","gaussian PDF",x,mean,sigma) ;
// Make extended PDF based on gauss. n will be the expected number of events
RooRealVar n("n","number of events",10000,0,20000) ;
RooExtendPdf egauss("egauss","extended gaussian PDF",gauss,n) ;
// Generate events from extended PDF
// The default number of events to generate is taken from gauss.expectedEvents()
// but can be overrided using a second argument
RooDataSet* data = egauss.generate(x) ;
// Fit PDF to dataset in extended mode (selected by fit option "e")
egauss.fitTo(*data,"mhe") ;
// Plot both on a frame ;
RooPlot* xframe = x.frame() ;
data->plotOn(xframe) ;
egauss->plotOn(xframe,Normalization(1.0,RooAbsReal::RelativeExpected)) ; // select intrinsic normalization
xframe->Draw() ;
// Make an extended gaussian PDF where the number of expected events
// is counted in a limited region of the dependent range
x.setRange("cut",-4,2) ;
RooRealVar mean2("mean2","mean of gaussian",-3) ;
RooRealVar sigma2("sigma2","width of gaussian",1) ;
RooGaussian gauss2("gauss2","gaussian PDF 2",x,mean2,sigma2) ;
RooRealVar n2("n2","number of events",10000,0,20000) ;
RooExtendPdf egauss2("egauss2","extended gaussian PDF w limited range",gauss2,n2,"cut") ;
egauss2.fitTo(*data,"mhe");
cout << "fitted number of events in data in range (-6,0) = " << n2.getVal() << endl ;
// Adding two extended PDFs gives an extended sum PDF
mean = 3.0 ; sigma = 2.0 ;
// Note that we omit coefficients when adding extended PDFS
RooAddPdf sumgauss("sumgauss","sum of two extended gauss PDFs",RooArgList(egauss,egauss2)) ;
sumgauss->plotOn(xframe,Normalization(1.0,RooAbsReal::RelativeExpected),LineColor(kRed)) ; // select intrinsic normalization
xframe->Draw() ;
// Note that in the plot sumgauss does not follow the normalization of the data
// because its expected number was intentionally chosen not to match the number of events in the data
// If no special 'cut normalizations' are needed (as done in egauss2), there is a shorthand
// way to construct an extended sumpdf:
RooAddPdf sumgauss2("sumgauss2","extended sum of two gaussian PDFs",
RooArgList(gauss,gauss2),RooArgList(n,n2)) ;
sumgauss2->plotOn(xframe,Normalization(1.0,RooAbsReal::RelativeExpected),LineColor(kGreen)) ; // select intrinsic normalization
xframe->Draw() ;
// Note that sumgauss2 looks different from sumgauss because for gauss2 the expected number
// of event parameter n2 now applies to the entire gauss2 area, whereas in egauss2 it was
// constructed to represent the number of events in the range (-4,-2). If we would use a separate
// parameter n3, set to 10000, to represent the number of events for gauss2 in sumgauss2, then
// sumgauss and sumgauss2 would be indentical.
}
Copyright © 2000-2005 University of California, Stanford University
|