63 ClassImp(RooStats::SPlot); ;
65 using namespace RooStats;
72 if(TestBit(kOwnData) && fSData)
93 SPlot::SPlot(
const char* name,
const char* title):
108 SPlot::SPlot(
const char* name,
const char* title,
const RooDataSet &data):
115 fSData = (RooDataSet*) &data;
121 SPlot::SPlot(
const SPlot &other):
124 RooArgList Args = (RooArgList) other.GetSWeightVars();
126 fSWeightVars.addClone(Args);
128 fSData = (RooDataSet*) other.GetSDataSet();
148 SPlot::SPlot(
const char* name,
const char* title, RooDataSet& data, RooAbsPdf* pdf,
149 const RooArgList &yieldsList,
const RooArgSet &projDeps,
150 bool includeWeights,
bool cloneData,
const char* newName,
151 const RooCmdArg& arg5,
const RooCmdArg& arg6,
const RooCmdArg& arg7,
const RooCmdArg& arg8):
155 fSData = (RooDataSet*) data.Clone(newName);
159 fSData = (RooDataSet*) &data;
162 for (
const auto arg : yieldsList) {
163 if (!dynamic_cast<const RooRealVar*>(arg)) {
164 coutE(InputArguments) <<
"SPlot::SPlot(" << GetName() <<
") input argument "
165 << arg->GetName() <<
" is not of type RooRealVar " << endl ;
166 throw string(Form(
"SPlot::SPlot(%s) input argument %s is not of type RooRealVar",GetName(),arg->GetName())) ;
174 this->AddSWeight(pdf, yieldsList, projDeps, includeWeights, arg5, arg6, arg7, arg8);
179 RooDataSet* SPlot::SetSData(RooDataSet* data)
182 fSData = (RooDataSet*) data;
191 RooDataSet* SPlot::GetSDataSet()
const
200 Double_t SPlot::GetSWeight(Int_t numEvent,
const char* sVariable)
const
202 if(numEvent > fSData->numEntries() )
204 coutE(InputArguments) <<
"Invalid Entry Number" << endl;
210 coutE(InputArguments) <<
"Invalid Entry Number" << endl;
214 Double_t totalYield = 0;
216 std::string varname(sVariable);
220 if(fSWeightVars.find(sVariable) )
222 RooArgSet Row(*fSData->get(numEvent));
223 totalYield += Row.getRealValue(sVariable);
228 if( fSWeightVars.find(varname.c_str()) )
231 RooArgSet Row(*fSData->get(numEvent));
232 totalYield += Row.getRealValue(varname.c_str() );
238 coutE(InputArguments) <<
"InputVariable not in list of sWeighted variables" << endl;
249 Double_t SPlot::GetSumOfEventSWeight(Int_t numEvent)
const
251 if(numEvent > fSData->numEntries() )
253 coutE(InputArguments) <<
"Invalid Entry Number" << endl;
259 coutE(InputArguments) <<
"Invalid Entry Number" << endl;
263 Int_t numSWeightVars = this->GetNumSWeightVars();
265 Double_t eventSWeight = 0;
267 RooArgSet Row(*fSData->get(numEvent));
269 for (Int_t i = 0; i < numSWeightVars; i++)
270 eventSWeight += Row.getRealValue(fSWeightVars.at(i)->GetName() );
280 Double_t SPlot::GetYieldFromSWeight(
const char* sVariable)
const
283 Double_t totalYield = 0;
285 std::string varname(sVariable);
289 if(fSWeightVars.find(sVariable) )
291 for(Int_t i=0; i < fSData->numEntries(); i++)
293 RooArgSet Row(*fSData->get(i));
294 totalYield += Row.getRealValue(sVariable);
300 if( fSWeightVars.find(varname.c_str()) )
302 for(Int_t i=0; i < fSData->numEntries(); i++)
304 RooArgSet Row(*fSData->get(i));
305 totalYield += Row.getRealValue(varname.c_str() );
312 coutE(InputArguments) <<
"InputVariable not in list of sWeighted variables" << endl;
321 RooArgList SPlot::GetSWeightVars()
const
324 RooArgList Args = fSWeightVars;
335 Int_t SPlot::GetNumSWeightVars()
const
337 RooArgList Args = fSWeightVars;
339 return Args.getSize();
368 void SPlot::AddSWeight( RooAbsPdf* pdf,
const RooArgList &yieldsTmp,
369 const RooArgSet &projDeps,
bool includeWeights,
370 const RooCmdArg& arg5,
const RooCmdArg& arg6,
const RooCmdArg& arg7,
const RooCmdArg& arg8)
373 RooFit::MsgLevel currentLevel = RooMsgService::instance().globalKillBelow();
377 RooArgList* constParameters = (RooArgList*)pdf->getParameters(fSData) ;
378 constParameters->remove(yieldsTmp, kTRUE, kTRUE);
383 std::vector<RooRealVar*> constVarHolder;
385 for(Int_t i = 0; i < constParameters->getSize(); i++)
387 RooRealVar* varTemp = (
dynamic_cast<RooRealVar*
>( constParameters->at(i) ) );
388 if(varTemp && varTemp->isConstant() == 0 )
390 varTemp->setConstant();
391 constVarHolder.push_back(varTemp);
398 pdf->fitTo(*fSData, RooFit::Extended(kTRUE), RooFit::SumW2Error(kTRUE), RooFit::PrintLevel(-1), RooFit::PrintEvalErrors(-1), arg5, arg6, arg7, arg8);
401 std::vector<double> yieldsHolder;
403 for(Int_t i = 0; i < yieldsTmp.getSize(); i++)
404 yieldsHolder.push_back( ((RooRealVar*) yieldsTmp.at(i))->getVal());
406 Int_t nspec = yieldsTmp.getSize();
407 RooArgList yields = *(RooArgList*)yieldsTmp.snapshot(kFALSE);
409 if(currentLevel <= RooFit::DEBUG)
411 coutI(InputArguments) <<
"Printing Yields" << endl;
417 RooArgSet vars(*fSData->get() );
418 vars.remove(projDeps, kTRUE, kTRUE);
424 pdf->attachDataSet(*fSData);
427 std::vector<RooRealVar*> yieldvars ;
428 RooArgSet* parameters = pdf->getParameters(fSData) ;
430 std::vector<Double_t> yieldvalues ;
431 for (Int_t k = 0; k < nspec; ++k)
433 RooRealVar* thisyield =
dynamic_cast<RooRealVar*
>(yields.at(k)) ;
435 RooRealVar* yieldinpdf =
dynamic_cast<RooRealVar*
>(parameters->find(thisyield->GetName() )) ;
438 coutI(InputArguments)<<
"yield in pdf: " << yieldinpdf->GetName() <<
" " << thisyield->getVal() << endl;
440 yieldvars.push_back(yieldinpdf) ;
441 yieldvalues.push_back(thisyield->getVal()) ;
446 Int_t numevents = fSData->numEntries() ;
448 std::vector<std::vector<Double_t> > pdfvalues(numevents,std::vector<Double_t>(nspec,0)) ;
452 for(Int_t m=0; m<nspec; ++m) yieldvars[m]->setVal(0) ;
461 RooArgSet * pdfvars = pdf->getVariables();
463 for (Int_t ievt = 0; ievt <numevents; ievt++)
471 RooStats::SetParameters(fSData->get(ievt), pdfvars);
475 for(Int_t k = 0; k < nspec; ++k)
478 if(yieldvars[k]->getMin() > 0)
480 coutW(InputArguments) <<
"Minimum Range for " << yieldvars[k]->GetName() <<
" must be 0. ";
481 coutW(InputArguments) <<
"Setting min range to 0" << std::endl;
482 yieldvars[k]->setMin(0);
485 if(yieldvars[k]->getMax() < 1)
487 coutW(InputArguments) <<
"Maximum Range for " << yieldvars[k]->GetName() <<
" must be 1. ";
488 coutW(InputArguments) <<
"Setting max range to 1" << std::endl;
489 yieldvars[k]->setMax(1);
493 yieldvars[k]->setVal( 1 ) ;
495 Double_t f_k = pdf->getVal(&vars) ;
496 pdfvalues[ievt][k] = f_k ;
497 if( !(f_k>1 || f_k<1) )
498 coutW(InputArguments) <<
"Strange pdf value: " << ievt <<
" " << k <<
" " << f_k << std::endl ;
499 yieldvars[k]->setVal( 0 ) ;
505 std::vector<Double_t> norm(nspec,0) ;
506 for (Int_t ievt = 0; ievt <numevents ; ievt++)
509 for(Int_t k=0; k<nspec; ++k) dnorm += yieldvalues[k] * pdfvalues[ievt][k] ;
510 for(Int_t j=0; j<nspec; ++j) norm[j] += pdfvalues[ievt][j]/dnorm ;
513 coutI(Contents) <<
"likelihood norms: " ;
515 for(Int_t k=0; k<nspec; ++k) coutI(Contents) << norm[k] <<
" " ;
516 coutI(Contents) << std::endl ;
519 TMatrixD covInv(nspec, nspec);
520 for (Int_t i = 0; i < nspec; i++)
for (Int_t j = 0; j < nspec; j++) covInv(i,j) = 0;
522 coutI(Contents) <<
"Calculating covariance matrix";
526 for (Int_t ievt = 0; ievt < numevents; ++ievt)
536 for(Int_t k = 0; k < nspec; ++k)
537 dsum += pdfvalues[ievt][k] * yieldvalues[k] ;
539 for(Int_t n=0; n<nspec; ++n)
540 for(Int_t j=0; j<nspec; ++j)
542 if(includeWeights == kTRUE)
543 covInv(n,j) += fSData->weight()*pdfvalues[ievt][n]*pdfvalues[ievt][j]/(dsum*dsum) ;
545 covInv(n,j) += pdfvalues[ievt][n]*pdfvalues[ievt][j]/(dsum*dsum) ;
555 if (covInv.Determinant() <=0)
557 coutE(Eval) <<
"SPlot Error: covariance matrix is singular; I can't invert it!" << std::endl;
562 TMatrixD covMatrix(TMatrixD::kInverted,covInv);
565 if(currentLevel <= RooFit::DEBUG)
567 coutI(Eval) <<
"Checking Likelihood normalization: " << std::endl;
568 coutI(Eval) <<
"Yield of specie Sum of Row in Matrix Norm" << std::endl;
569 for(Int_t k=0; k<nspec; ++k)
571 Double_t covnorm(0) ;
572 for(Int_t m=0; m<nspec; ++m) covnorm += covInv[k][m]*yieldvalues[m] ;
574 for(Int_t m = 0; m < nspec; ++m) sumrow += covMatrix[k][m] ;
575 coutI(Eval) << yieldvalues[k] <<
" " << sumrow <<
" " << covnorm << endl ;
580 coutI(Eval) <<
"Calculating sWeight" << std::endl;
581 std::vector<RooRealVar*> sweightvec ;
582 std::vector<RooRealVar*> pdfvec ;
583 RooArgSet sweightset ;
588 fSWeightVars.Clear();
590 for(Int_t k=0; k<nspec; ++k)
592 std::string wname = std::string(yieldvars[k]->GetName()) +
"_sw";
593 RooRealVar* var =
new RooRealVar(wname.c_str(),wname.c_str(),0) ;
594 sweightvec.push_back( var) ;
595 sweightset.add(*var) ;
596 fSWeightVars.add(*var);
598 wname =
"L_" + std::string(yieldvars[k]->GetName());
599 var =
new RooRealVar(wname.c_str(),wname.c_str(),0) ;
600 pdfvec.push_back( var) ;
601 sweightset.add(*var) ;
607 RooDataSet* sWeightData =
new RooDataSet(
"dataset",
"dataset with sWeights", sweightset);
609 for(Int_t ievt = 0; ievt < numevents; ++ievt)
616 for(Int_t k = 0; k < nspec; ++k) dsum += pdfvalues[ievt][k] * yieldvalues[k] ;
618 for(Int_t n=0; n<nspec; ++n)
621 for(Int_t j=0; j<nspec; ++j) nsum += covMatrix(n,j) * pdfvalues[ievt][j] ;
629 if(includeWeights == kTRUE) sweightvec[n]->setVal(fSData->weight() * nsum/dsum) ;
630 else sweightvec[n]->setVal( nsum/dsum) ;
632 pdfvec[n]->setVal( pdfvalues[ievt][n] ) ;
634 if( !(fabs(nsum/dsum)>=0 ) )
636 coutE(Contents) <<
"error: " << nsum/dsum << endl ;
641 sWeightData->add(sweightset) ;
647 fSData->merge(sWeightData);
653 for(Int_t i = 0; i < yieldsTmp.getSize(); i++)
654 ((RooRealVar*) yieldsTmp.at(i))->setVal(yieldsHolder.at(i));
658 for(Int_t i=0; i < (Int_t) constVarHolder.size(); i++)
659 constVarHolder.at(i)->setConstant(kFALSE);