Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
SPlot.h
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer 21/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#ifndef RooStats_SPlot
13#define RooStats_SPlot
14
15class RooAbsReal;
16class RooAbsPdf;
17class RooFitResult;
18class RooRealVar;
19class RooSimultaneous;
20
21
22#include "RooMsgService.h"
23
24#include "RooFitResult.h"
25#include "RooRealVar.h"
26#include "RooHist.h"
27#include "RooPlot.h"
28#include "RooDataSet.h"
29
30namespace RooStats{
31
32 class SPlot: public TNamed {
33
34 public:
35
36 ~SPlot() override;
37 SPlot();
38 SPlot(const SPlot &other);
39 SPlot(const char* name, const char* title);
40 SPlot(const char* name, const char* title, const RooDataSet &data);
41 SPlot(const char* name, const char* title,RooDataSet& data, RooAbsPdf* pdf,
42 const RooArgList &yieldsList,const RooArgSet &projDeps=RooArgSet(),
43 bool useWeights=true, bool copyDataSet = false, const char* newName = "",
44 const RooCmdArg& fitToarg5={},
45 const RooCmdArg& fitToarg6={},
46 const RooCmdArg& fitToarg7={},
47 const RooCmdArg& fitToarg8={});
48
50
51 RooDataSet* GetSDataSet() const;
52
54
56
57 void AddSWeight(RooAbsPdf* pdf, const RooArgList &yieldsTmp,
58 const RooArgSet &projDeps=RooArgSet(), bool includeWeights=true,
59 const RooCmdArg& fitToarg5={},
60 const RooCmdArg& fitToarg6={},
61 const RooCmdArg& fitToarg7={},
62 const RooCmdArg& fitToarg8={});
63
64 double GetSumOfEventSWeight(Int_t numEvent) const;
65
66 double GetYieldFromSWeight(const char* sVariable) const;
67
68 double GetSWeight(Int_t numEvent, const char* sVariable) const;
69
70
71
72 protected:
73
74 enum {
75 kOwnData = BIT(20)
76 };
77
79
80 // RooListProxy fSWeightVars;
81
82 RooDataSet *fSData = nullptr;
83
84 ClassDefOverride(SPlot,1) // Class used for making sPlots
85
86
87 };
88
89}
90#endif
#define BIT(n)
Definition Rtypes.h:85
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
char name[80]
Definition TGX11.cxx:110
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
Container class to hold unbinned data.
Definition RooDataSet.h:57
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
A class to calculate "sWeights" used to create an "sPlot".
Definition SPlot.h:32
double GetSWeight(Int_t numEvent, const char *sVariable) const
Retrieve an s weight.
Definition SPlot.cxx:236
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.
Definition SPlot.cxx:404
RooDataSet * fSData
Definition SPlot.h:82
RooArgList GetSWeightVars() const
Return a RooArgList containing all parameters that have s weights.
Definition SPlot.cxx:357
double GetSumOfEventSWeight(Int_t numEvent) const
Sum the SWeights for a particular event.
Definition SPlot.cxx:285
SPlot()
Default constructor.
Definition SPlot.cxx:130
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:371
~SPlot() override
Definition SPlot.cxx:120
RooDataSet * SetSData(RooDataSet *data)
Set dataset (if not passed in constructor).
Definition SPlot.cxx:215
RooDataSet * GetSDataSet() const
Retrieve s-weighted data.
Definition SPlot.cxx:227
double GetYieldFromSWeight(const char *sVariable) const
Sum the SWeights for a particular species over all events.
Definition SPlot.cxx:316
RooArgList fSWeightVars
Definition SPlot.h:78
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
Namespace for the RooStats classes.
Definition Asimov.h:19