Logo ROOT  
Reference Guide
SPlot.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer 28/07/2008
3
4/*************************************************************************
5 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/*****************************************************************************
13 * Project: RooStats
14 * Package: RooFit/RooStats
15 *
16 * Authors:
17 * Original code from M. Pivk as part of MLFit package from BaBar.
18 * Modifications:
19 * Giacinto Piacquadio, Maurizio Pierini: modifications for new RooFit version
20 * George H. Lewis, Kyle Cranmer: generalized for weighted events
21 *
22 * Porting to RooStats (with permission) by Kyle Cranmer, July 2008
23 *
24 *****************************************************************************/
25
26
27/** \class RooStats::SPlot
28 \ingroup Roostats
29
30 A class to calculate "sWeights" used to create an "sPlot".
31 An sPlot can reweight a dataset to show different components (e.g. signal / background),
32 but it doesn't use cuts, and therefore doesn't have to sort events into signal/background (or other) categories.
33 Instead of *assigning* a category to each event in the dataset, all events are *weighted*.
34 To compute the weights, a PDF with different components is analysed, and the weights are added
35 to the dataset. When plotting the dataset with the weights of the signal or
36 background components, the data looks like "signal", but all events in the dataset are used.
37
38 The result is similar to a likelihood projection plot, but without cuts.
39
40 \note SPlot needs to fit the pdf to the data once, so make sure that all relevant fit arguments such as
41 the fit range are passed in the constructor.
42
43 The code is based on
44 ``SPlot: A statistical tool to unfold data distributions,''
45 Nucl. Instrum. Meth. A 555, 356 (2005) [arXiv:physics/0402083].
46
47 ### Creating an SPlot
48 To use this class, you first must have a pdf that includes
49 yield parameters for (possibly several) different species, for example a signal and background
50 yield. Those yields must be of type RooRealVar / RooLinearVar (or anything that derives from
51 RooAbsRealLValue). This is necessary because
52 RooStats needs to be able to set the yields to 0 and 1 to probe the PDF. After
53 constructing the s weights, the yields will be restored to their original values.
54
55 To create an instance of the SPlot, supply a data set, the pdf to analyse,
56 and a list which parameters of the pdf are yields. The SPlot will calculate SWeights, and
57 include these as columns in the RooDataSet. The dataset will have two additional columns
58 for every yield with name "`<varname>`":
59 - `L_<varname>` is the the likelihood for each event, *i.e.*, the pdf evaluated for the given value of the variable "varname".
60 - `<varname>_sw` is the value of the sWeight for the variable "varname" for each event.
61
62 In SPlot::SPlot(), one can choose whether columns should be added to an existing dataset or whether a copy of the dataset
63 should be created.
64
65 ### Plotting s-weighted data
66 After computing the s weights, create a new dataset that uses the s weights of the variable of interest for weighting.
67 If the yield parameter for signal was e.g. "signalYield", the dataset can be constructed as follows:
68 ~~~{.cpp}
69 RooDataSet data_signal("<name>", "<title>", <dataWithSWeights>, <variables>, 0, "signalYield_sw");
70 ~~~
71
72 A complete tutorial with an extensive model is rs301_splot.C
73
74 #### Using ratios as yield parameters
75 As mentioned, RooStats needs to be able to modify the yield parameters. That means that they have to be a RooRealVar
76 of a RooLinearVar. This allows using ratio parameters as in the following example:
77 ~~~{.cpp}
78 RooRealVar x("x", "observable", 0, 0, 20);
79 RooRealVar m("m", "mean", 5., -10, 10);
80 RooRealVar s("s", "sigma", 2., 0, 10);
81 RooGaussian gaus("gaus", "gaus", x, m, s);
82
83 RooRealVar a("a", "exp", -0.2, -10., 0.);
84 RooExponential ex("ex", "ex", x, a);
85
86 RooRealVar common("common", "common scale", 3., 0, 10);
87 RooRealVar r1("r1", "ratio of signal events", 0.3, 0, 10);
88 RooRealVar r2("r2", "ratio of background events", 0.5, 0, 10);
89 RooLinearVar c1("c1", "c1", r1, common, RooFit::RooConst(0.));
90 RooLinearVar c2("c2", "c2", r2, common, RooFit::RooConst(0.));
91
92 RooAddPdf sum("sum", "sum", RooArgSet(gaus, ex), RooArgSet(c1, c2));
93 auto data = sum.generate(x, 1000);
94
95 RooStats::SPlot splot("splot", "splot", *data, &sum, RooArgSet(c1, c2));
96 ~~~
97*/
98
99#include <vector>
100#include <map>
101
102#include "RooStats/SPlot.h"
103#include "RooAbsPdf.h"
104#include "RooDataSet.h"
105#include "RooRealVar.h"
106#include "RooGlobalFunc.h"
107#include "TTree.h"
109
110
111#include "TMatrixD.h"
112
113
115
116using namespace RooStats;
117using namespace std;
118
119////////////////////////////////////////////////////////////////////////////////
120
121SPlot::~SPlot()
122{
123 if(TestBit(kOwnData) && fSData)
124 delete fSData;
125
126}
127
128////////////////////////////////////////////////////////////////////////////////
129/// Default constructor
130
132 TNamed()
133{
134 RooArgList Args;
135
136 fSWeightVars = Args;
137
138 fSData = NULL;
139
140}
141
142////////////////////////////////////////////////////////////////////////////////
143
144SPlot::SPlot(const char* name, const char* title):
145 TNamed(name, title)
146{
147 RooArgList Args;
148
149 fSWeightVars = Args;
150
151 fSData = NULL;
152
153}
154
155////////////////////////////////////////////////////////////////////////////////
156///Constructor from a RooDataSet
157///No sWeighted variables are present
158
159SPlot::SPlot(const char* name, const char* title, const RooDataSet &data):
160 TNamed(name, title)
161{
162 RooArgList Args;
163
164 fSWeightVars = Args;
165
166 fSData = (RooDataSet*) &data;
167}
168
169////////////////////////////////////////////////////////////////////////////////
170/// Copy Constructor from another SPlot
171
172SPlot::SPlot(const SPlot &other):
173 TNamed(other)
174{
175 RooArgList Args = (RooArgList) other.GetSWeightVars();
176
178
179 fSData = (RooDataSet*) other.GetSDataSet();
180
181}
182
183////////////////////////////////////////////////////////////////////////////////
184///Construct a new SPlot instance, calculate sWeights, and include them
185///in the RooDataSet held by this instance.
186///
187/// The constructor automatically calls AddSWeight() to add s weights to the dataset.
188/// These can be retrieved later using GetSWeight() or GetSDataSet().
189///\param[in] name Name of the instance.
190///\param[in] title Title of the instance.
191///\param[in] data Dataset to fit to.
192///\param[in] pdf PDF to compute s weights for.
193///\param[in] yieldsList List of parameters in `pdf` that are yields. These must be RooRealVar or RooLinearVar, since RooStats will need to modify their values.
194///\param[in] projDeps Don't normalise over these parameters when calculating the sWeights. Will be passed on to AddSWeight().
195///\param[in] useWeights Include weights of the input data in calculation of s weights.
196///\param[in] cloneData Make a clone of the incoming data before adding weights.
197///\param[in] newName New name for the data.
198///\param[in] argX Additional arguments for the fitting step in AddSWeight().
199SPlot::SPlot(const char* name, const char* title, RooDataSet& data, RooAbsPdf* pdf,
200 const RooArgList &yieldsList, const RooArgSet &projDeps,
201 bool useWeights, bool cloneData, const char* newName,
202 const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8):
203 TNamed(name, title)
204{
205 if(cloneData == 1) {
206 fSData = (RooDataSet*) data.Clone(newName);
208 }
209 else
210 fSData = (RooDataSet*) &data;
211
212 // Add check that yieldsList contains all RooRealVar / RooAbsRealLValue
213 for (const auto arg : yieldsList) {
214 if (!dynamic_cast<const RooAbsRealLValue*>(arg)) {
215 coutE(InputArguments) << "SPlot::SPlot(" << GetName() << ") input argument "
216 << arg->GetName() << " is not of type RooRealVar (or RooLinearVar)."
217 << "\nRooStats must be able to set it to 0 and to 1 to probe the PDF." << endl ;
218 throw std::invalid_argument(Form("SPlot::SPlot(%s) input argument %s is not of type RooRealVar/RooLinearVar",GetName(),arg->GetName())) ;
219 }
220 }
221
222 //Construct a new SPlot class,
223 //calculate sWeights, and include them
224 //in the RooDataSet of this class.
225
226 this->AddSWeight(pdf, yieldsList, projDeps, useWeights, arg5, arg6, arg7, arg8);
227}
228
229////////////////////////////////////////////////////////////////////////////////
230/// Set dataset (if not passed in constructor).
232{
233 if(data) {
234 fSData = (RooDataSet*) data;
235 return fSData;
236 } else
237 return NULL;
238}
239
240////////////////////////////////////////////////////////////////////////////////
241/// Retrieve s-weighted data.
242/// It does **not** automatically call AddSWeight(). This needs to be done manually.
244{
245 return fSData;
246}
247
248////////////////////////////////////////////////////////////////////////////////
249/// Retrieve an s weight.
250/// \param[in] numEvent Event number to retrieve s weight for.
251/// \param[in] sVariable The yield parameter to retrieve the s weight for.
252Double_t SPlot::GetSWeight(Int_t numEvent, const char* sVariable) const
253{
254 if(numEvent > fSData->numEntries() )
255 {
256 coutE(InputArguments) << "Invalid Entry Number" << endl;
257 return -1;
258 }
259
260 if(numEvent < 0)
261 {
262 coutE(InputArguments) << "Invalid Entry Number" << endl;
263 return -1;
264 }
265
266 Double_t totalYield = 0;
267
268 std::string varname(sVariable);
269 varname += "_sw";
270
271
272 if(fSWeightVars.find(sVariable) )
273 {
274 RooArgSet Row(*fSData->get(numEvent));
275 totalYield += Row.getRealValue(sVariable);
276
277 return totalYield;
278 }
279
280 if( fSWeightVars.find(varname.c_str()) )
281 {
282
283 RooArgSet Row(*fSData->get(numEvent));
284 totalYield += Row.getRealValue(varname.c_str() );
285
286 return totalYield;
287 }
288
289 else
290 coutE(InputArguments) << "InputVariable not in list of sWeighted variables" << endl;
291
292 return -1;
293}
294
295
296////////////////////////////////////////////////////////////////////////////////
297/// Sum the SWeights for a particular event.
298/// This sum should equal the total weight of that event.
299/// This method is intended to be used as a check.
300
302{
303 if(numEvent > fSData->numEntries() )
304 {
305 coutE(InputArguments) << "Invalid Entry Number" << endl;
306 return -1;
307 }
308
309 if(numEvent < 0)
310 {
311 coutE(InputArguments) << "Invalid Entry Number" << endl;
312 return -1;
313 }
314
315 Int_t numSWeightVars = this->GetNumSWeightVars();
316
317 Double_t eventSWeight = 0;
318
319 RooArgSet Row(*fSData->get(numEvent));
320
321 for (Int_t i = 0; i < numSWeightVars; i++)
322 eventSWeight += Row.getRealValue(fSWeightVars.at(i)->GetName() );
323
324 return eventSWeight;
325}
326
327////////////////////////////////////////////////////////////////////////////////
328/// Sum the SWeights for a particular species over all events.
329/// This should equal the total (weighted) yield of that species.
330/// This method is intended as a check.
331
332Double_t SPlot::GetYieldFromSWeight(const char* sVariable) const
333{
334
335 Double_t totalYield = 0;
336
337 std::string varname(sVariable);
338 varname += "_sw";
339
340
341 if(fSWeightVars.find(sVariable) )
342 {
343 for(Int_t i=0; i < fSData->numEntries(); i++)
344 {
345 RooArgSet Row(*fSData->get(i));
346 totalYield += Row.getRealValue(sVariable);
347 }
348
349 return totalYield;
350 }
351
352 if( fSWeightVars.find(varname.c_str()) )
353 {
354 for(Int_t i=0; i < fSData->numEntries(); i++)
355 {
356 RooArgSet Row(*fSData->get(i));
357 totalYield += Row.getRealValue(varname.c_str() );
358 }
359
360 return totalYield;
361 }
362
363 else
364 coutE(InputArguments) << "InputVariable not in list of sWeighted variables" << endl;
365
366 return -1;
367}
368
369
370////////////////////////////////////////////////////////////////////////////////
371/// Return a RooArgList containing all paramters that have s weights.
372
374{
375
377
378 return Args;
379
380}
381
382////////////////////////////////////////////////////////////////////////////////
383/// Return the number of SWeights
384/// In other words, return the number of
385/// species that we are trying to extract.
386
388{
390
391 return Args.getSize();
392}
393
394////////////////////////////////////////////////////////////////////////////////
395/// Method which adds the sWeights to the dataset.
396///
397/// The SPlot will contain two new variables for each yield parameter:
398/// - `L_<varname>` is the the likelihood for each event, *i.e.*, the pdf evaluated for the a given value of the variable "varname".
399/// - `<varname>_sw` is the value of the sWeight for the variable "varname" for each event.
400///
401/// Find Parameters in the PDF to be considered fixed when calculating the SWeights
402/// and be sure to NOT include the yields in that list.
403///
404/// After fixing non-yield parameters, this function will start a fit by calling
405/// ```
406/// pdf->fitTo(*fSData, RooFit::Extended(kTRUE), RooFit::SumW2Error(kTRUE), RooFit::PrintLevel(-1), RooFit::PrintEvalErrors(-1)).
407/// ```
408/// One can pass additional arguments to `fitTo`, such as `RooFit::Range("fitrange")`, as `arg5`, `arg6`, `arg7`, `arg8`.
409///
410/// \note A `RooFit::Range` may be necessary to get expected results if you initially fit in a range
411/// and/or called `pdf->fixCoefRange("fitrange")` on `pdf`.
412/// Pass `arg5`, `arg6`, `arg7`, `arg8` AT YOUR OWN RISK.
413///
414/// \param[in] pdf PDF to fit to data to compute s weights.
415/// \param[in] yieldsTmp Yields to use to compute s weights.
416/// \param[in] projDeps These will not be normalized over when calculating the sWeights,
417/// and will be considered parameters, not observables.
418/// \param[in] includeWeights Include weights of the input data in calculation of s weights.
419/// \param[in] argX Optional additional arguments for the fitting step.
420void SPlot::AddSWeight( RooAbsPdf* pdf, const RooArgList &yieldsTmp,
421 const RooArgSet &projDeps, bool includeWeights,
422 const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
423{
424
425 // Find Parameters in the PDF to be considered fixed when calculating the SWeights
426 // and be sure to NOT include the yields in that list
427 RooArgList* constParameters = (RooArgList*)pdf->getParameters(fSData) ;
428 for (unsigned int i=0; i < constParameters->size(); ++i) {
429 // Need a counting loop since collection is being modified
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, kTRUE, kTRUE);
433 --i;
434 }
435 }
436
437
438 // Set these parameters constant and store them so they can later
439 // be set to not constant
440 std::vector<RooAbsRealLValue*> constVarHolder;
441
442 for(Int_t i = 0; i < constParameters->getSize(); i++)
443 {
444 RooAbsRealLValue* varTemp = static_cast<RooAbsRealLValue*>( constParameters->at(i) );
445 if(varTemp && varTemp->isConstant() == 0 )
446 {
447 varTemp->setConstant();
448 constVarHolder.push_back(varTemp);
449 }
450 }
451
452 // Fit yields to the data with all other variables held constant
453 // This is necessary because SPlot assumes the yields minimise -Log(likelihood)
455
456 // Hold the value of the fitted yields
457 std::vector<double> yieldsHolder;
458
459 for(Int_t i = 0; i < yieldsTmp.getSize(); i++)
460 yieldsHolder.push_back(static_cast<RooAbsReal*>(yieldsTmp.at(i))->getVal());
461
462 const Int_t nspec = yieldsTmp.getSize();
463 RooArgList yields = *(RooArgList*)yieldsTmp.snapshot(kFALSE);
464
466 coutI(InputArguments) << "Printing Yields" << endl;
467 yields.Print();
468 }
469
470 // The list of variables to normalize over when calculating PDF values.
471
472 RooArgSet vars(*fSData->get() );
473 vars.remove(projDeps, kTRUE, kTRUE);
474
475 // Attach data set
476
477 // const_cast<RooAbsPdf*>(pdf)->attachDataSet(*fSData);
478
479 pdf->attachDataSet(*fSData);
480
481 // first calculate the pdf values for all species and all events
482 std::vector<RooAbsRealLValue*> yieldvars ;
483 RooArgSet pdfServers;
484 pdf->treeNodeServerList(&pdfServers);
485
486 std::vector<Double_t> yieldvalues ;
487 for (Int_t k = 0; k < nspec; ++k) {
488 auto thisyield = static_cast<const RooAbsReal*>(yields.at(k)) ;
489 auto yieldinpdf = static_cast<RooAbsRealLValue*>(pdfServers.find(thisyield->GetName()));
490 assert(pdf->dependsOn(*yieldinpdf));
491
492 if (yieldinpdf) {
493 coutI(InputArguments)<< "yield in pdf: " << yieldinpdf->GetName() << " " << thisyield->getVal() << endl;
494
495 yieldvars.push_back(yieldinpdf) ;
496 yieldvalues.push_back(thisyield->getVal()) ;
497 }
498 }
499
500 Int_t numevents = fSData->numEntries() ;
501
502
503
504
505 // set all yield to zero
506 for(Int_t m=0; m<nspec; ++m) {
507 auto theVar = static_cast<RooAbsRealLValue*>(yieldvars[m]);
508 theVar->setVal(0) ;
509
510 //Check that range of yields is at least (0,1), and fix otherwise
511 if (theVar->getMin() > 0) {
512 coutE(InputArguments) << "Yield variables need to have a range that includes at least [0, 1]. Minimum for "
513 << theVar->GetName() << " is " << theVar->getMin() << std::endl;
514 if (RooRealVar* realVar = dynamic_cast<RooRealVar*>(theVar)) {
515 coutE(InputArguments) << "Setting min range to 0" << std::endl;
516 realVar->setMin(0);
517 } else {
518 throw std::invalid_argument(std::string("Yield variable ") + theVar->GetName() + " must have a range that includes 0.");
519 }
520 }
521
522 if (theVar->getMax() < 1) {
523 coutW(InputArguments) << "Yield variables need to have a range that includes at least [0, 1]. Maximum for "
524 << theVar->GetName() << " is " << theVar->getMax() << std::endl;
525 if (RooRealVar* realVar = dynamic_cast<RooRealVar*>(theVar)) {
526 coutE(InputArguments) << "Setting max range to 1" << std::endl;
527 realVar->setMax(1);
528 } else {
529 throw std::invalid_argument(std::string("Yield variable ") + theVar->GetName() + " must have a range that includes 1.");
530 }
531 }
532 }
533
534
535 // For every event and for every species,
536 // calculate the value of the component pdf for that specie
537 // by setting the yield of that specie to 1
538 // and all others to 0. Evaluate the pdf for each event
539 // and store the values.
540
541 RooArgSet * pdfvars = pdf->getVariables();
542 std::vector<std::vector<Double_t> > pdfvalues(numevents,std::vector<Double_t>(nspec,0)) ;
543
544 for (Int_t ievt = 0; ievt <numevents; ievt++)
545 {
546 //WVE: FIX THIS PART, EVALUATION PROGRESS!!
547
548 RooStats::SetParameters(fSData->get(ievt), pdfvars);
549
550 for(Int_t k = 0; k < nspec; ++k) {
551 auto theVar = static_cast<RooAbsRealLValue*>(yieldvars[k]);
552
553 // set this yield to 1
554 theVar->setVal( 1 ) ;
555 // evaluate the pdf
556 Double_t f_k = pdf->getVal(&vars) ;
557 pdfvalues[ievt][k] = f_k ;
558 if( !(f_k>1 || f_k<1) )
559 coutW(InputArguments) << "Strange pdf value: " << ievt << " " << k << " " << f_k << std::endl ;
560 theVar->setVal( 0 ) ;
561 }
562 }
563 delete pdfvars;
564
565 // check that the likelihood normalization is fine
566 std::vector<Double_t> norm(nspec,0) ;
567 for (Int_t ievt = 0; ievt <numevents ; ievt++)
568 {
569 Double_t dnorm(0) ;
570 for(Int_t k=0; k<nspec; ++k) dnorm += yieldvalues[k] * pdfvalues[ievt][k] ;
571 for(Int_t j=0; j<nspec; ++j) norm[j] += pdfvalues[ievt][j]/dnorm ;
572 }
573
574 coutI(Contents) << "likelihood norms: " ;
575
576 for(Int_t k=0; k<nspec; ++k) coutI(Contents) << norm[k] << " " ;
577 coutI(Contents) << std::endl ;
578
579 // Make a TMatrixD to hold the covariance matrix.
580 TMatrixD covInv(nspec, nspec);
581 for (Int_t i = 0; i < nspec; i++) for (Int_t j = 0; j < nspec; j++) covInv(i,j) = 0;
582
583 coutI(Contents) << "Calculating covariance matrix";
584
585
586 // Calculate the inverse covariance matrix, using weights
587 for (Int_t ievt = 0; ievt < numevents; ++ievt)
588 {
589
590 fSData->get(ievt) ;
591
592 // Calculate contribution to the inverse of the covariance
593 // matrix. See BAD 509 V2 eqn. 15
594
595 // Sum for the denominator
596 Double_t dsum(0);
597 for(Int_t k = 0; k < nspec; ++k)
598 dsum += pdfvalues[ievt][k] * yieldvalues[k] ;
599
600 for(Int_t n=0; n<nspec; ++n)
601 for(Int_t j=0; j<nspec; ++j)
602 {
603 if(includeWeights)
604 covInv(n,j) += fSData->weight()*pdfvalues[ievt][n]*pdfvalues[ievt][j]/(dsum*dsum) ;
605 else
606 covInv(n,j) += pdfvalues[ievt][n]*pdfvalues[ievt][j]/(dsum*dsum) ;
607 }
608
609 //ADDED WEIGHT ABOVE
610
611 }
612
613 // Covariance inverse should now be computed!
614
615 // Invert to get the covariance matrix
616 if (covInv.Determinant() <=0)
617 {
618 coutE(Eval) << "SPlot Error: covariance matrix is singular; I can't invert it!" << std::endl;
619 covInv.Print();
620 return;
621 }
622
623 TMatrixD covMatrix(TMatrixD::kInverted,covInv);
624
625 //check cov normalization
626 if (RooMsgService::instance().isActive(this, RooFit::Eval, RooFit::DEBUG)) {
627 coutI(Eval) << "Checking Likelihood normalization: " << std::endl;
628 coutI(Eval) << "Yield of specie Sum of Row in Matrix Norm" << std::endl;
629 for(Int_t k=0; k<nspec; ++k)
630 {
631 Double_t covnorm(0) ;
632 for(Int_t m=0; m<nspec; ++m) covnorm += covInv[k][m]*yieldvalues[m] ;
633 Double_t sumrow(0) ;
634 for(Int_t m = 0; m < nspec; ++m) sumrow += covMatrix[k][m] ;
635 coutI(Eval) << yieldvalues[k] << " " << sumrow << " " << covnorm << endl ;
636 }
637 }
638
639 // calculate for each event the sWeight (BAD 509 V2 eq. 21)
640 coutI(Eval) << "Calculating sWeight" << std::endl;
641 std::vector<RooRealVar*> sweightvec ;
642 std::vector<RooRealVar*> pdfvec ;
643 RooArgSet sweightset ;
644
645 // Create and label the variables
646 // used to store the SWeights
647
649
650 for(Int_t k=0; k<nspec; ++k)
651 {
652 std::string wname = std::string(yieldvars[k]->GetName()) + "_sw";
653 RooRealVar* var = new RooRealVar(wname.c_str(),wname.c_str(),0) ;
654 sweightvec.push_back( var) ;
655 sweightset.add(*var) ;
656 fSWeightVars.add(*var);
657
658 wname = "L_" + std::string(yieldvars[k]->GetName());
659 var = new RooRealVar(wname.c_str(),wname.c_str(),0) ;
660 pdfvec.push_back( var) ;
661 sweightset.add(*var) ;
662 }
663
664 // Create and fill a RooDataSet
665 // with the SWeights
666
667 RooDataSet* sWeightData = new RooDataSet("dataset", "dataset with sWeights", sweightset);
668
669 for(Int_t ievt = 0; ievt < numevents; ++ievt)
670 {
671
672 fSData->get(ievt) ;
673
674 // sum for denominator
675 Double_t dsum(0);
676 for(Int_t k = 0; k < nspec; ++k) dsum += pdfvalues[ievt][k] * yieldvalues[k] ;
677 // covariance weighted pdf for each specief
678 for(Int_t n=0; n<nspec; ++n)
679 {
680 Double_t nsum(0) ;
681 for(Int_t j=0; j<nspec; ++j) nsum += covMatrix(n,j) * pdfvalues[ievt][j] ;
682
683
684 //Add the sWeights here!!
685 //Include weights,
686 //ie events weights are absorbed into sWeight
687
688
689 if(includeWeights) sweightvec[n]->setVal(fSData->weight() * nsum/dsum) ;
690 else sweightvec[n]->setVal( nsum/dsum) ;
691
692 pdfvec[n]->setVal( pdfvalues[ievt][n] ) ;
693
694 if( !(fabs(nsum/dsum)>=0 ) )
695 {
696 coutE(Contents) << "error: " << nsum/dsum << endl ;
697 return;
698 }
699 }
700
701 sWeightData->add(sweightset) ;
702 }
703
704
705 // Add the SWeights to the original data set
706
707 fSData->merge(sWeightData);
708
709 delete sWeightData;
710
711 //Restore yield values
712
713 for(Int_t i = 0; i < yieldsTmp.getSize(); i++)
714 static_cast<RooAbsRealLValue*>(yieldsTmp.at(i))->setVal(yieldsHolder.at(i));
715
716 //Make any variables that were forced to constant no longer constant
717
718 for(Int_t i=0; i < (Int_t) constVarHolder.size(); i++)
719 constVarHolder.at(i)->setConstant(kFALSE);
720
721 return;
722
723}
#define coutI(a)
Definition: RooMsgService.h:30
#define coutW(a)
Definition: RooMsgService.h:32
#define coutE(a)
Definition: RooMsgService.h:33
int Int_t
Definition: RtypesCore.h:43
const Bool_t kFALSE
Definition: RtypesCore.h:90
const Bool_t kTRUE
Definition: RtypesCore.h:89
#define ClassImp(name)
Definition: Rtypes.h:361
char name[80]
Definition: TGX11.cxx:109
char * Form(const char *fmt,...)
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:73
Bool_t dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0, Bool_t valueOnly=kFALSE) const
Test whether we depend on (ie, are served by) any object in the specified collection.
Definition: RooAbsArg.cxx:730
void treeNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=0, Bool_t doBranch=kTRUE, Bool_t doLeaf=kTRUE, Bool_t valueOnly=kFALSE, Bool_t recurseNonDerived=kFALSE) const
Fill supplied list with nodes of the arg tree, following all server links, starting with ourself as t...
Definition: RooAbsArg.cxx:502
RooArgSet * getParameters(const RooAbsData *data, Bool_t stripDisconnected=kTRUE) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
Definition: RooAbsArg.cxx:544
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:1909
Bool_t isConstant() const
Check if the "Constant" attribute is set.
Definition: RooAbsArg.h:342
void attachDataSet(const RooAbsData &set)
Replace server nodes with names matching the dataset variable names with those data set variables,...
Definition: RooAbsArg.cxx:1467
Int_t getSize() const
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents.
virtual RooAbsArg * addClone(const RooAbsArg &var, Bool_t silent=kFALSE)
Add a clone of the specified argument to list.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
const_iterator end() const
Double_t getRealValue(const char *name, Double_t defVal=0, Bool_t verbose=kFALSE) const
Get value of a RooAbsReal stored in set with given name.
Storage_t::size_type size() const
const_iterator begin() const
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:306
virtual RooFitResult * fitTo(RooAbsData &data, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Fit PDF to given dataset.
Definition: RooAbsPdf.cxx:1254
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
virtual void setVal(Double_t value)=0
void setConstant(Bool_t value=kTRUE)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:90
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition: RooArgList.h:74
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition: RooCmdArg.h:28
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:33
virtual TObject * Clone(const char *newname="") const override
Make a clone of an object using the Streamer facility.
Definition: RooDataSet.h:68
virtual const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
Bool_t merge(RooDataSet *data1, RooDataSet *data2=0, RooDataSet *data3=0, RooDataSet *data4=0, RooDataSet *data5=0, RooDataSet *data6=0)
virtual void add(const RooArgSet &row, Double_t weight=1.0, Double_t weightError=0) override
Add a data point, with its coordinates specified in the 'data' argset, to the data set.
virtual Double_t weight() const override
Return event weight of current event.
Definition: RooDataSet.cxx:973
static RooMsgService & instance()
Return reference to singleton instance.
Bool_t isActive(const RooAbsArg *self, RooFit::MsgTopic facility, RooFit::MsgLevel level)
Check if logging is active for given object/topic/RooFit::MsgLevel combination.
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
A class to calculate "sWeights" used to create an "sPlot".
Definition: SPlot.h:32
void AddSWeight(RooAbsPdf *pdf, const RooArgList &yieldsTmp, const RooArgSet &projDeps=RooArgSet(), bool includeWeights=kTRUE, const RooCmdArg &fitToarg5=RooCmdArg::none(), const RooCmdArg &fitToarg6=RooCmdArg::none(), const RooCmdArg &fitToarg7=RooCmdArg::none(), const RooCmdArg &fitToarg8=RooCmdArg::none())
Method which adds the sWeights to the dataset.
Definition: SPlot.cxx:420
Double_t GetYieldFromSWeight(const char *sVariable) const
Sum the SWeights for a particular species over all events.
Definition: SPlot.cxx:332
Double_t GetSWeight(Int_t numEvent, const char *sVariable) const
Retrieve an s weight.
Definition: SPlot.cxx:252
RooDataSet * fSData
Definition: SPlot.h:82
RooArgList GetSWeightVars() const
Return a RooArgList containing all paramters that have s weights.
Definition: SPlot.cxx:373
SPlot()
Default constructor.
Definition: SPlot.cxx:131
Int_t GetNumSWeightVars() const
Return the number of SWeights In other words, return the number of species that we are trying to extr...
Definition: SPlot.cxx:387
RooDataSet * SetSData(RooDataSet *data)
Set dataset (if not passed in constructor).
Definition: SPlot.cxx:231
RooDataSet * GetSDataSet() const
Retrieve s-weighted data.
Definition: SPlot.cxx:243
RooArgList fSWeightVars
Definition: SPlot.h:78
Double_t GetSumOfEventSWeight(Int_t numEvent) const
Sum the SWeights for a particular event.
Definition: SPlot.cxx:301
void Print(Option_t *name="") const
Print the matrix as a table of elements.
virtual Double_t Determinant() const
Return the matrix determinant.
Definition: TMatrixT.cxx:1364
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual void Clear(Option_t *="")
Definition: TObject.h:115
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:187
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
RooCmdArg SumW2Error(Bool_t flag)
RooCmdArg Extended(Bool_t flag=kTRUE)
RooCmdArg PrintEvalErrors(Int_t numErrors)
RooCmdArg PrintLevel(Int_t code)
const Int_t n
Definition: legend1.C:16
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
@ InputArguments
Definition: RooGlobalFunc.h:68
Namespace for the RooStats classes.
Definition: Asimov.h:19
void SetParameters(const RooArgSet *desiredVals, RooArgSet *paramsToChange)
Definition: RooStatsUtils.h:65
auto * m
Definition: textangle.C:8