Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
RooMultiPdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * RW, Ruddick William UC Colorado wor@slac.stanford.edu *
7 * *
8 * Copyright (c) 2000-2005, Regents of the University of California *
9 * and Stanford University. All rights reserved. *
10 * *
11 * Redistribution and use in source and binary forms, *
12 * with or without modification, are permitted according to the terms *
13 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
14 *****************************************************************************/
15/** \class RooMultiPdf
16 \ingroup Roofit
17
18 The class RooMultiPdf allows for the creation of a RooMultiPdf object,
19 which can switch between previously set probability density functions (PDF)
20 by setting their associated indices.*/
21
22#include <RooMultiPdf.h>
24#include <RooConstVar.h>
25
26// Constructing a RooMultiPdf
27// parameter name : The name of the RooMultiPdf object
28// parameter title : Display title in plots
29// parameter _x : variable used to select which PDF that selects the active PDF.
30// parameter _c : A list of the pdfs.The index of each PDF in the list should match the values in _x
31RooMultiPdf::RooMultiPdf(const char *name, const char *title, RooCategory &_x, const RooArgList &_c)
32 : RooAbsPdf(name, title), // call of constructor base class RooAbsPdf passing it name and title
33 c("_pdfs", "The list of pdfs", this),
34 corr("_corrs", "The list of correction factors", this),
35 x("_index", "the pdf index", this, _x)
36
37// parameter corr : Holds correction factors - number of free parameters
38// in each PDF held in the RooMultiPdf object
39{
40 int count = 0;
41
42 c.add(_c);
43 for (RooAbsArg *pdf : c) {
44 // This is done by the user BUT is there a way to do it at construction?
45 _x.defineType(("_pdf" + std::to_string(count)).c_str(), count);
46 std::unique_ptr<RooArgSet> variables(pdf->getVariables());
47 std::unique_ptr<RooAbsCollection> nonConstVariables(variables->selectByAttrib("Constant", false));
48 // Isn't there a better way to hold on to these values?
49 std::string corrName = std::string{"const"} + pdf->GetName();
50 corr.addOwned(std::make_unique<RooConstVar>(corrName.c_str(), "", nonConstVariables->size()));
51 count++;
52 }
54}
55
56// Here new RooMultiPdf copy is created that references the same components as the original.
57// Copies c, corr, x
58RooMultiPdf::RooMultiPdf(const RooMultiPdf &other, const char *name)
59 : RooAbsPdf(other, name), c("_pdfs", this, other.c), corr("_corrs", this, other.corr), x("_index", this, other.x)
60{
61 fIndex = other.fIndex;
63 cFactor = other.cFactor; // correction to 2*NLL by default is -> 2*0.5 per param
64}
65
66// evaluate() and getLogVal() define how the value and log-value of the RooMultiPdf are computed at a given point.
67// RooMultiPdf must have both of these so it RooFit can treat it like a PDF
69{
70 double val = getCurrentPdf()->getVal(c.nset());
71 _oldIndex = x;
72 return val;
73}
74
76{
77 double logval = getCurrentPdf()->getLogVal(nset);
78 _oldIndex = x;
79 return logval;
80}
81
82void RooMultiPdf::getParametersHook(const RooArgSet *nset, RooArgSet *list, bool stripDisconnected) const
83{
84 if (!stripDisconnected)
85 return;
86
87 list->removeAll();
88 getCurrentPdf()->getParameters(nset, *list, stripDisconnected);
89 list->add(*x);
90}
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
char name[80]
Definition TGX11.cxx:148
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...
RooAbsArg()
Default constructor.
RooAbsPdf()
Default constructor.
virtual double getLogVal(const RooArgSet *set=nullptr) const
Return the log of the current value with given normalization An error message is printed if the argum...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
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:24
Object to represent discrete states.
Definition RooCategory.h:28
bool defineType(const std::string &label)
Define a state with given name.
RooCategoryProxy x
Definition RooMultiPdf.h:40
RooListProxy c
Definition RooMultiPdf.h:38
void getParametersHook(const RooArgSet *nset, RooArgSet *list, bool stripDisconnected) const override
RooAbsPdf * getCurrentPdf() const
Definition RooMultiPdf.h:22
RooListProxy corr
Definition RooMultiPdf.h:39
Double_t evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Double_t getLogVal(const RooArgSet *set=nullptr) const override
Return the log of the current value with given normalization An error message is printed if the argum...
double cFactor
Definition RooMultiPdf.h:48
Int_t _oldIndex
Definition RooMultiPdf.h:44