185 bool useWeights,
bool cloneData,
const char* newName,
197 for (
const auto arg : yieldsList) {
199 coutE(InputArguments) <<
"SPlot::SPlot(" <<
GetName() <<
") input argument "
200 << arg->GetName() <<
" is not of type RooRealVar (or RooLinearVar)."
201 <<
"\nRooStats must be able to set it to 0 and to 1 to probe the PDF." << endl ;
202 throw std::invalid_argument(
Form(
"SPlot::SPlot(%s) input argument %s is not of type RooRealVar/RooLinearVar",
GetName(),arg->GetName())) ;
210 this->
AddSWeight(pdf, yieldsList, projDeps, useWeights, arg5, arg6, arg7, arg8);
240 coutE(InputArguments) <<
"Invalid Entry Number" << endl;
246 coutE(InputArguments) <<
"Invalid Entry Number" << endl;
250 double totalYield = 0;
252 std::string varname(sVariable);
274 coutE(InputArguments) <<
"InputVariable not in list of sWeighted variables" << endl;
289 coutE(InputArguments) <<
"Invalid Entry Number" << endl;
295 coutE(InputArguments) <<
"Invalid Entry Number" << endl;
301 double eventSWeight = 0;
305 for (
Int_t i = 0; i < numSWeightVars; i++)
319 double totalYield = 0;
321 std::string varname(sVariable);
348 coutE(InputArguments) <<
"InputVariable not in list of sWeighted variables" << endl;
405 const RooArgSet &projDeps,
bool includeWeights,
412 for (
unsigned int i=0; i < constParameters->size(); ++i) {
414 auto& par = *(*constParameters)[i];
415 if (std::any_of(yieldsTmp.
begin(), yieldsTmp.
end(), [&](
const RooAbsArg* yield){ return yield->dependsOn(par); })) {
416 constParameters->remove(par,
true,
true);
424 std::vector<RooAbsRealLValue*> constVarHolder;
426 for(std::size_t i = 0; i < constParameters->size(); i++)
432 constVarHolder.push_back(varTemp);
442 vars.
remove(projDeps,
true,
true);
445 std::vector<double> yieldsHolder;
447 yieldsHolder.reserve(yieldsTmp.
size());
448 for(std::size_t i = 0; i < yieldsTmp.
size(); i++) {
456 coutI(InputArguments) <<
"Printing Yields" << endl;
462 std::vector<RooAbsRealLValue*> yieldvars ;
466 std::vector<double> yieldvalues ;
467 for (
Int_t k = 0; k < nspec; ++k) {
468 auto thisyield =
static_cast<const RooAbsReal*
>(yields.
at(k)) ;
473 coutI(InputArguments)<<
"yield in pdf: " << yieldinpdf->GetName() <<
" " << thisyield->getVal(&vars) << endl;
475 yieldvars.push_back(yieldinpdf) ;
476 yieldvalues.push_back(thisyield->getVal(&vars)) ;
491 if (theVar->getMin() > 0) {
492 coutE(InputArguments) <<
"Yield variables need to have a range that includes at least [0, 1]. Minimum for "
493 << theVar->GetName() <<
" is " << theVar->getMin() << std::endl;
495 coutE(InputArguments) <<
"Setting min range to 0" << std::endl;
498 throw std::invalid_argument(std::string(
"Yield variable ") + theVar->GetName() +
" must have a range that includes 0.");
502 if (theVar->getMax() < 1) {
503 coutW(InputArguments) <<
"Yield variables need to have a range that includes at least [0, 1]. Maximum for "
504 << theVar->GetName() <<
" is " << theVar->getMax() << std::endl;
506 coutE(InputArguments) <<
"Setting max range to 1" << std::endl;
509 throw std::invalid_argument(std::string(
"Yield variable ") + theVar->GetName() +
" must have a range that includes 1.");
521 std::unique_ptr<RooArgSet> pdfvars{pdf->
getVariables()};
522 std::vector<std::vector<double> > pdfvalues(numevents,std::vector<double>(nspec,0)) ;
524 for (
Int_t ievt = 0; ievt <numevents; ievt++)
530 for(
Int_t k = 0; k < nspec; ++k) {
536 double f_k = pdf->
getVal(&vars) ;
537 pdfvalues[ievt][k] = f_k ;
538 if( !(f_k>1 || f_k<1) )
539 coutW(InputArguments) <<
"Strange pdf value: " << ievt <<
" " << k <<
" " << f_k << std::endl ;
540 theVar->setVal( 0 ) ;
545 std::vector<double> norm(nspec,0) ;
546 for (
Int_t ievt = 0; ievt <numevents ; ievt++)
549 for(
Int_t k=0; k<nspec; ++k) dnorm += yieldvalues[k] * pdfvalues[ievt][k] ;
550 for(
Int_t j=0; j<nspec; ++j) norm[j] += pdfvalues[ievt][j]/dnorm ;
553 coutI(Contents) <<
"likelihood norms: " ;
555 for(
Int_t k=0; k<nspec; ++k)
coutI(Contents) << norm[k] <<
" " ;
556 coutI(Contents) << std::endl ;
560 for (
Int_t i = 0; i < nspec; i++)
for (
Int_t j = 0; j < nspec; j++) covInv(i,j) = 0;
562 coutI(Contents) <<
"Calculating covariance matrix";
566 for (
Int_t ievt = 0; ievt < numevents; ++ievt)
576 for(
Int_t k = 0; k < nspec; ++k)
577 dsum += pdfvalues[ievt][k] * yieldvalues[k] ;
580 for(
Int_t j=0; j<nspec; ++j)
582 if (includeWeights) {
583 covInv(
n, j) +=
fSData->
weight() * pdfvalues[ievt][
n] * pdfvalues[ievt][j] / (dsum * dsum);
585 covInv(
n, j) += pdfvalues[ievt][
n] * pdfvalues[ievt][j] / (dsum * dsum);
599 coutE(Eval) <<
"SPlot Error: covariance matrix is singular; I can't invert it!" << std::endl;
608 coutI(Eval) <<
"Checking Likelihood normalization: " << std::endl;
609 coutI(Eval) <<
"Yield of specie Sum of Row in Matrix Norm" << std::endl;
610 for(
Int_t k=0; k<nspec; ++k)
613 for(
Int_t m=0;
m<nspec; ++
m) covnorm += covInv[k][
m]*yieldvalues[
m] ;
615 for(
Int_t m = 0;
m < nspec; ++
m) sumrow += covMatrix[k][
m] ;
616 coutI(Eval) << yieldvalues[k] <<
" " << sumrow <<
" " << covnorm << endl ;
621 coutI(Eval) <<
"Calculating sWeight" << std::endl;
622 std::vector<RooRealVar*> sweightvec ;
623 std::vector<RooRealVar*> pdfvec ;
631 for(
Int_t k=0; k<nspec; ++k)
633 std::string wname = std::string(yieldvars[k]->
GetName()) +
"_sw";
635 sweightvec.push_back( var) ;
636 sweightset.
add(*var) ;
639 wname =
"L_" + std::string(yieldvars[k]->
GetName());
640 var =
new RooRealVar(wname.c_str(),wname.c_str(),0) ;
641 pdfvec.push_back( var) ;
642 sweightset.
add(*var) ;
650 for(
Int_t ievt = 0; ievt < numevents; ++ievt)
657 for(
Int_t k = 0; k < nspec; ++k) dsum += pdfvalues[ievt][k] * yieldvalues[k] ;
662 for(
Int_t j=0; j<nspec; ++j) nsum += covMatrix(
n,j) * pdfvalues[ievt][j] ;
670 if(includeWeights) sweightvec[
n]->setVal(
fSData->
weight() * nsum/dsum) ;
671 else sweightvec[
n]->setVal( nsum/dsum) ;
673 pdfvec[
n]->setVal( pdfvalues[ievt][
n] ) ;
675 if( !(std::abs(nsum/dsum)>=0 ) )
677 coutE(Contents) <<
"error: " << nsum/dsum << endl ;
682 sWeightData->
add(sweightset) ;
694 for(std::size_t i = 0; i < yieldsTmp.
size(); i++)
699 for(
Int_t i=0; i < (
Int_t) constVarHolder.size(); i++)
700 constVarHolder.at(i)->setConstant(
false);
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Common abstract base class for objects that represent a value and a "shape" in RooFit.
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
bool isConstant() const
Check if the "Constant" attribute is set.
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
void treeNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=nullptr, bool doBranch=true, bool doLeaf=true, bool valueOnly=false, bool recurseNonDerived=false) const
Fill supplied list with nodes of the arg tree, following all server links, starting with ourself as t...
double getRealValue(const char *name, double defVal=0.0, bool verbose=false) const
Get value of a RooAbsReal stored in set with given name.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
RooAbsCollection * snapshot(bool deepCopy=true) const
Take a snap shot of current collection contents.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
const_iterator end() const
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
Storage_t::size_type size() const
virtual RooAbsArg * addClone(const RooAbsArg &var, bool silent=false)
Add a clone of the specified argument to list.
const_iterator begin() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract interface for all probability density functions.
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
void setConstant(bool value=true)
virtual void setVal(double value)=0
Set the current value of the object. Needs to be overridden by implementations.
Abstract base class for objects that represent a real value and implements functionality common to al...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
RooArgList is a container object that can hold multiple RooAbsArg objects.
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Named container for two doubles, two integers two object points and three string pointers that can be...
Container class to hold unbinned data.
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
bool merge(RooDataSet *data1, RooDataSet *data2=nullptr, RooDataSet *data3=nullptr, RooDataSet *data4=nullptr, RooDataSet *data5=nullptr, RooDataSet *data6=nullptr)
void add(const RooArgSet &row, double weight, double weightError)
Add one ore more rows of data.
double weight() const override
Return event weight of current event.
static RooMsgService & instance()
Return reference to singleton instance.
bool isActive(T self, RooFit::MsgTopic topic, RooFit::MsgLevel level)
Check if logging is active for given object/topic/RooFit::MsgLevel combination.
Variable that can be changed from the outside.
A class to calculate "sWeights" used to create an "sPlot".
double GetSWeight(Int_t numEvent, const char *sVariable) const
Retrieve an s weight.
void AddSWeight(RooAbsPdf *pdf, const RooArgList &yieldsTmp, const RooArgSet &projDeps=RooArgSet(), bool includeWeights=true, const RooCmdArg &fitToarg5={}, const RooCmdArg &fitToarg6={}, const RooCmdArg &fitToarg7={}, const RooCmdArg &fitToarg8={})
Method which adds the sWeights to the dataset.
RooArgList GetSWeightVars() const
Return a RooArgList containing all parameters that have s weights.
double GetSumOfEventSWeight(Int_t numEvent) const
Sum the SWeights for a particular event.
SPlot()
Default constructor.
Int_t GetNumSWeightVars() const
Return the number of SWeights In other words, return the number of species that we are trying to extr...
RooDataSet * SetSData(RooDataSet *data)
Set dataset (if not passed in constructor).
RooDataSet * GetSDataSet() const
Retrieve s-weighted data.
double GetYieldFromSWeight(const char *sVariable) const
Sum the SWeights for a particular species over all events.
void Print(Option_t *name="") const override
Print the matrix as a table of elements.
Double_t Determinant() const override
Return the matrix determinant.
The TNamed class is the base class for all named ROOT classes.
const char * GetName() const override
Returns name of object.
virtual void Clear(Option_t *="")
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
RooCmdArg SumW2Error(bool flag)
RooCmdArg PrintEvalErrors(Int_t numErrors)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Extended(bool flag=true)
Namespace for the RooStats classes.