121 if(TestBit(kOwnData) && fSData)
133 fSWeightVars.
assign(Args);
141 fSWeightVars.
assign(Args);
153 fSWeightVars.
assign(Args);
184 bool useWeights,
bool cloneData,
const char* newName,
189 fSData =
static_cast<RooDataSet*
>(data.Clone(newName));
196 for (
const auto arg : yieldsList) {
198 coutE(InputArguments) <<
"SPlot::SPlot(" << GetName() <<
") input argument "
199 << arg->GetName() <<
" is not of type RooRealVar (or RooLinearVar)."
200 <<
"\nRooStats must be able to set it to 0 and to 1 to probe the PDF." << std::endl ;
201 throw std::invalid_argument(
Form(
"SPlot::SPlot(%s) input argument %s is not of type RooRealVar/RooLinearVar",GetName(),arg->GetName())) ;
209 this->AddSWeight(pdf, yieldsList, projDeps, useWeights, arg5, arg6, arg7, arg8);
237 if(numEvent > fSData->numEntries() )
239 coutE(InputArguments) <<
"Invalid Entry Number" << std::endl;
245 coutE(InputArguments) <<
"Invalid Entry Number" << std::endl;
249 double totalYield = 0;
251 std::string varname(sVariable);
255 if(fSWeightVars.find(sVariable) )
258 totalYield += Row.getRealValue(sVariable);
263 if( fSWeightVars.find(varname.c_str()) )
267 totalYield += Row.getRealValue(varname.c_str() );
273 coutE(InputArguments) <<
"InputVariable not in list of sWeighted variables" << std::endl;
286 if(numEvent > fSData->numEntries() )
288 coutE(InputArguments) <<
"Invalid Entry Number" << std::endl;
294 coutE(InputArguments) <<
"Invalid Entry Number" << std::endl;
298 Int_t numSWeightVars = this->GetNumSWeightVars();
300 double eventSWeight = 0;
304 for (
Int_t i = 0; i < numSWeightVars; i++)
305 eventSWeight += Row.getRealValue(fSWeightVars.at(i)->GetName() );
318 double totalYield = 0;
320 std::string varname(sVariable);
324 if(fSWeightVars.find(sVariable) )
326 for(
Int_t i=0; i < fSData->numEntries(); i++)
329 totalYield += Row.getRealValue(sVariable);
335 if( fSWeightVars.find(varname.c_str()) )
337 for(
Int_t i=0; i < fSData->numEntries(); i++)
340 totalYield += Row.getRealValue(varname.c_str() );
347 coutE(InputArguments) <<
"InputVariable not in list of sWeighted variables" << std::endl;
410 std::unique_ptr<RooArgSet> discrimVars{pdf->
getObservables(fSData)};
411 for (
auto const *yield : yieldsTmp) {
412 if (discrimVars->find(yield->GetName()) || yield->dependsOn(*discrimVars)) {
413 coutE(InputArguments) <<
"SPlot::AddSWeight(" << GetName() <<
") the yield variable \"" << yield->GetName()
414 <<
"\" is also a discriminating variable used in the fit, "
415 <<
"or depends on one. The sPlot formalism requires the yield parameters to be "
416 <<
"independent of the discriminating observables. The resulting sWeights would "
417 <<
"be invalid. Please remove this variable from either the fit observables or "
418 <<
"from the list of yields." << std::endl;
419 throw std::invalid_argument(
Form(
"SPlot::AddSWeight(%s) yield variable \"%s\" is also used "
420 "as a discriminating variable in the fit.",
421 GetName(), yield->GetName()));
427 std::unique_ptr<RooArgSet> constParameters{pdf->
getParameters(fSData)};
428 for (
unsigned int i=0; i < constParameters->size(); ++i) {
430 auto& par = *(*constParameters)[i];
431 if (std::any_of(yieldsTmp.begin(), yieldsTmp.end(), [&](
const RooAbsArg* yield){ return yield->dependsOn(par); })) {
432 constParameters->remove(par,
true,
true);
440 std::vector<RooAbsRealLValue*> constVarHolder;
442 for(std::size_t i = 0; i < constParameters->size(); i++)
448 constVarHolder.push_back(varTemp);
458 vars.remove(projDeps,
true,
true);
461 std::vector<double> yieldsHolder;
463 yieldsHolder.reserve(yieldsTmp.size());
464 for(std::size_t i = 0; i < yieldsTmp.size(); i++) {
465 yieldsHolder.push_back(
static_cast<RooAbsReal*
>(yieldsTmp.at(i))->
getVal(&vars));
468 const Int_t nspec = yieldsTmp.size();
472 coutI(InputArguments) <<
"Printing Yields" << std::endl;
478 std::vector<RooAbsRealLValue*> yieldvars ;
482 std::vector<double> yieldvalues ;
483 for (
Int_t k = 0; k < nspec; ++k) {
484 auto thisyield =
static_cast<const RooAbsReal*
>(yields.
at(k)) ;
489 coutI(InputArguments)<<
"yield in pdf: " << yieldinpdf->GetName() <<
" " << thisyield->getVal(&vars) << std::endl;
491 yieldvars.push_back(yieldinpdf) ;
492 yieldvalues.push_back(thisyield->getVal(&vars)) ;
496 Int_t numevents = fSData->numEntries() ;
507 if (theVar->getMin() > 0) {
508 coutE(InputArguments) <<
"Yield variables need to have a range that includes at least [0, 1]. Minimum for "
509 << theVar->GetName() <<
" is " << theVar->getMin() << std::endl;
511 coutE(InputArguments) <<
"Setting min range to 0" << std::endl;
514 throw std::invalid_argument(std::string(
"Yield variable ") + theVar->GetName() +
" must have a range that includes 0.");
518 if (theVar->getMax() < 1) {
519 coutW(InputArguments) <<
"Yield variables need to have a range that includes at least [0, 1]. Maximum for "
520 << theVar->GetName() <<
" is " << theVar->getMax() << std::endl;
522 coutE(InputArguments) <<
"Setting max range to 1" << std::endl;
525 throw std::invalid_argument(std::string(
"Yield variable ") + theVar->GetName() +
" must have a range that includes 1.");
537 std::unique_ptr<RooArgSet> pdfvars{pdf->
getVariables()};
538 std::vector<std::vector<double> > pdfvalues(numevents,std::vector<double>(nspec,0)) ;
540 for (
Int_t ievt = 0; ievt <numevents; ievt++)
544 pdfvars->assign(*fSData->get(ievt));
546 for(
Int_t k = 0; k < nspec; ++k) {
552 double f_k = pdf->
getVal(&vars) ;
553 pdfvalues[ievt][k] = f_k ;
554 if( !(f_k>1 || f_k<1) )
555 coutW(InputArguments) <<
"Strange pdf value: " << ievt <<
" " << k <<
" " << f_k << std::endl ;
556 theVar->setVal( 0 ) ;
561 std::vector<double> norm(nspec,0) ;
562 for (
Int_t ievt = 0; ievt <numevents ; ievt++)
565 for(
Int_t k=0; k<nspec; ++k) dnorm += yieldvalues[k] * pdfvalues[ievt][k] ;
566 for(
Int_t j=0; j<nspec; ++j) norm[j] += pdfvalues[ievt][j]/dnorm ;
569 coutI(Contents) <<
"likelihood norms: " ;
571 for(
Int_t k=0; k<nspec; ++k)
coutI(Contents) << norm[k] <<
" " ;
572 coutI(Contents) << std::endl ;
576 for (
Int_t i = 0; i < nspec; i++)
for (
Int_t j = 0; j < nspec; j++) covInv(i,j) = 0;
578 coutI(Contents) <<
"Calculating covariance matrix";
582 for (
Int_t ievt = 0; ievt < numevents; ++ievt)
592 for(
Int_t k = 0; k < nspec; ++k)
593 dsum += pdfvalues[ievt][k] * yieldvalues[k] ;
596 for(
Int_t j=0; j<nspec; ++j)
598 if (includeWeights) {
599 covInv(
n, j) += fSData->weight() * pdfvalues[ievt][
n] * pdfvalues[ievt][j] / (dsum * dsum);
601 covInv(
n, j) += pdfvalues[ievt][
n] * pdfvalues[ievt][j] / (dsum * dsum);
613 if (covInv.Determinant() <=0)
615 coutE(Eval) <<
"SPlot Error: covariance matrix is singular; I can't invert it!" << std::endl;
624 coutI(Eval) <<
"Checking Likelihood normalization: " << std::endl;
625 coutI(Eval) <<
"Yield of specie Sum of Row in Matrix Norm" << std::endl;
626 for(
Int_t k=0; k<nspec; ++k)
629 for(
Int_t m=0;
m<nspec; ++
m) covnorm += covInv[k][
m]*yieldvalues[
m] ;
631 for(
Int_t m = 0;
m < nspec; ++
m) sumrow += covMatrix[k][
m] ;
632 coutI(Eval) << yieldvalues[k] <<
" " << sumrow <<
" " << covnorm << std::endl ;
637 coutI(Eval) <<
"Calculating sWeight" << std::endl;
638 std::vector<RooRealVar*> sweightvec ;
639 std::vector<RooRealVar*> pdfvec ;
645 fSWeightVars.
Clear();
647 for(
Int_t k=0; k<nspec; ++k)
649 std::string wname = std::string(yieldvars[k]->GetName()) +
"_sw";
651 sweightvec.push_back( var) ;
652 sweightset.
add(*var) ;
653 fSWeightVars.add(*var);
655 wname =
"L_" + std::string(yieldvars[k]->GetName());
656 var =
new RooRealVar(wname.c_str(),wname.c_str(),0) ;
657 pdfvec.push_back( var) ;
658 sweightset.
add(*var) ;
666 for(
Int_t ievt = 0; ievt < numevents; ++ievt)
673 for(
Int_t k = 0; k < nspec; ++k) dsum += pdfvalues[ievt][k] * yieldvalues[k] ;
678 for(
Int_t j=0; j<nspec; ++j) nsum += covMatrix(
n,j) * pdfvalues[ievt][j] ;
686 if(includeWeights) sweightvec[
n]->setVal(fSData->weight() * nsum/dsum) ;
687 else sweightvec[
n]->setVal( nsum/dsum) ;
689 pdfvec[
n]->setVal( pdfvalues[ievt][
n] ) ;
691 if( !(std::abs(nsum/dsum)>=0 ) )
693 coutE(Contents) <<
"error: " << nsum/dsum << std::endl ;
698 sWeightData->
add(sweightset) ;
704 fSData->merge(sWeightData);
710 for(std::size_t i = 0; i < yieldsTmp.size(); i++)
715 for(
Int_t i=0; i < (
Int_t) constVarHolder.size(); i++)
716 constVarHolder.at(i)->setConstant(
false);
int Int_t
Signed integer 4 bytes (int).
TMatrixT< Double_t > TMatrixD
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 > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
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...
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
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.
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.
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.
void add(const RooArgSet &row, double weight, double weightError)
Add one ore more rows of data.
static RooMsgService & instance()
Return reference to singleton instance.
Variable that can be changed from the outside.
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={})
Int_t GetNumSWeightVars() const
RooArgList GetSWeightVars() const
double GetSWeight(Int_t numEvent, const char *sVariable) const
double GetYieldFromSWeight(const char *sVariable) const
double GetSumOfEventSWeight(Int_t numEvent) const
RooDataSet * GetSDataSet() const
RooDataSet * SetSData(RooDataSet *data)
The TNamed class is the base class for all named ROOT classes.
virtual void Clear(Option_t *="")
RooCmdArg SumW2Error(bool flag)
RooCmdArg PrintEvalErrors(Int_t numErrors)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Extended(bool flag=true)
Namespace for the RooStats classes.