Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooResolutionModel.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * File: $Id: RooResolutionModel.h,v 1.26 2007/05/14 18:37:46 wouter Exp $
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16#ifndef ROO_RESOLUTION_MODEL
17#define ROO_RESOLUTION_MODEL
18
19#include "RooAbsPdf.h"
20#include "RooTemplateProxy.h"
21#include "RooRealVar.h"
22#include "RooFormulaVar.h"
23
25
27public:
28
29 // Constructors, assignment etc
30 inline RooResolutionModel() = default;
31 RooResolutionModel(const char *name, const char *title, RooAbsRealLValue& x) ;
32 RooResolutionModel(const RooResolutionModel& other, const char* name=nullptr);
33 TObject* clone(const char* newname) const override = 0;
34 ~RooResolutionModel() override;
35
37 const RooDataSet*, const RooArgSet*,
38 bool) const { return nullptr; }
39
40 double getValV(const RooArgSet* nset=nullptr) const override ;
41
42 // If used as regular PDF, it also has to be normalized. If this resolution
43 // model is used in a convolution, return unnormalized value regardless of
44 // specified normalization set.
45 bool selfNormalized() const override { return isConvolved() ; }
46
48 /// Return the convolution variable of the resolution model.
49 RooAbsRealLValue& convVar() const {return *x;}
50 const RooRealVar& basisConvVar() const ;
51
52 inline bool isBasisSupported(const char* name) const { return basisCode(name)?true:false ; }
53 virtual Int_t basisCode(const char* name) const = 0 ;
54
55 virtual void normLeafServerList(RooArgSet& list) const ;
56 double getNorm(const RooArgSet* nset=nullptr) const override ;
57
58 inline const RooFormulaVar& basis() const { return _basis?*_basis:*identity() ; }
59 bool isConvolved() const { return _basis ? true : false ; }
60
61 void printMultiline(std::ostream& os, Int_t content, bool verbose=false, TString indent="") const override ;
62
63 static RooFormulaVar* identity() ;
64
65 virtual void changeBasis(RooFormulaVar* basis) ;
66
67protected:
68
69 friend class RooConvGenContext ;
70 friend class RooAddModel ;
71 RooTemplateProxy<RooAbsRealLValue> x; ///< Dependent/convolution variable
72
73 bool redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive) override ;
74
75 friend class RooAbsAnaConvPdf ;
76
77 Int_t _basisCode ; ///< Identifier code for selected basis function
78 RooFormulaVar* _basis = nullptr; ///< Basis function convolved with this resolution model
79 bool _ownBasis ; ///< Flag indicating ownership of _basis
80
81 ClassDefOverride(RooResolutionModel, 2) // Abstract Resolution Model
82};
83
84#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition TGX11.cxx:110
Base class for PDFs that represent a physics model that can be analytically convolved with a resoluti...
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
Abstract container object that can hold multiple RooAbsArg objects.
Abstract base class for generator contexts of RooAbsPdf objects.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
RooAddModel is an efficient implementation of a sum of PDFs of the form.
Definition RooAddModel.h:27
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Efficient implementation of the generator context specific for RooAbsAnaConvPdf objects.
Container class to hold unbinned data.
Definition RooDataSet.h:34
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Variable that can be changed from the outside.
Definition RooRealVar.h:37
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
~RooResolutionModel() override
Destructor.
TObject * clone(const char *newname) const override=0
double getValV(const RooArgSet *nset=nullptr) const override
Modified version of RooAbsPdf::getValF().
virtual void changeBasis(RooFormulaVar *basis)
Change the basis function we convolute with.
virtual Int_t basisCode(const char *name) const =0
virtual RooResolutionModel * convolution(RooFormulaVar *basis, RooAbsArg *owner) const
Instantiate a clone of this resolution model representing a convolution with given basis function.
double getNorm(const RooArgSet *nset=nullptr) const override
Return the integral of this PDF over all elements of 'nset'.
static RooFormulaVar * identity()
Return identity formula pointer.
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive) override
Forward redirectServers call to our basis function, which is not connected to either resolution model...
const RooRealVar & basisConvVar() const
Return the convolution variable of the selection basis function.
bool _ownBasis
Flag indicating ownership of _basis.
bool isBasisSupported(const char *name) const
bool selfNormalized() const override
Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
virtual void normLeafServerList(RooArgSet &list) const
Floating point error checking and tracing for given float value.
virtual RooAbsGenContext * modelGenContext(const RooAbsAnaConvPdf &, const RooArgSet &, const RooDataSet *, const RooArgSet *, bool) const
void printMultiline(std::ostream &os, Int_t content, bool verbose=false, TString indent="") const override
Print info about this object to the specified stream.
Int_t _basisCode
Identifier code for selected basis function.
RooAbsRealLValue & convVar() const
Return the convolution variable of the resolution model.
RooResolutionModel()=default
RooFormulaVar * _basis
Basis function convolved with this resolution model.
const RooFormulaVar & basis() const
RooTemplateProxy< RooAbsRealLValue > x
Dependent/convolution variable.
Mother of all ROOT objects.
Definition TObject.h:41
Basic string class.
Definition TString.h:139