Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
PiecewiseInterpolation.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * File: $Id: PiecewiseInterpolation.h,v 1.3 2007/05/11 09:11:30 verkerke 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_PIECEWISEINTERPOLATION
17#define ROO_PIECEWISEINTERPOLATION
18
19#include "RooAbsReal.h"
20#include "RooRealProxy.h"
21#include "RooListProxy.h"
22
23#include "RooObjCacheManager.h"
24#include <vector>
25#include <list>
26
27class RooRealVar;
28class RooArgList;
29
31public:
32
34 PiecewiseInterpolation(const char *name, const char *title, const RooAbsReal &nominal, const RooArgList &lowSet,
35 const RooArgList &highSet, const RooArgList &paramSet);
36 ~PiecewiseInterpolation() override ;
37
38 PiecewiseInterpolation(const PiecewiseInterpolation& other, const char *name = nullptr);
39 TObject* clone(const char* newname) const override { return new PiecewiseInterpolation(*this, newname); }
40
41 /// Return pointer to the nominal hist function.
42 const RooAbsReal* nominalHist() const {
43 return &_nominal.arg();
44 }
45
46 // virtual double defaultErrorLevel() const ;
47
48 // void printMetaArgs(std::ostream& os) const ;
49
50 const RooArgList& lowList() const { return _lowSet ; }
51 const RooArgList& highList() const { return _highSet ; }
52 const RooArgList& paramList() const { return _paramSet ; }
53 const std::vector<int>& interpolationCodes() const { return _interpCode; }
54
55 //virtual bool forceAnalyticalInt(const RooAbsArg&) const { return true ; }
56 bool setBinIntegrator(RooArgSet& allVars) ;
57
58 Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet,const char* rangeName=nullptr) const override ;
59 double analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=nullptr) const override ;
60
61 void setPositiveDefinite(bool flag=true){_positiveDefinite=flag;}
62 bool positiveDefinite() const {return _positiveDefinite;}
63
64 void setInterpCode(RooAbsReal& param, int code, bool silent=false);
65 void setAllInterpCodes(int code);
67
68 std::list<double>* binBoundaries(RooAbsRealLValue& /*obs*/, double /*xlo*/, double /*xhi*/) const override ;
69 std::list<double>* plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const override ;
70 bool isBinnedDistribution(const RooArgSet& obs) const override ;
71
72 void translate(RooFit::Detail::CodeSquashContext &ctx) const override;
73
74protected:
75
77 public:
78 CacheElem() {} ;
81 ret.add(_lowIntList);
82 ret.add(_highIntList);
83 return ret ;
84 }
88 // will want std::vector<RooRealVar*> for low and high also
89 } ;
90 mutable RooObjCacheManager _normIntMgr ; ///<! The integration cache manager
91
92 RooRealProxy _nominal; ///< The nominal value
93 RooArgList _ownedList ; ///< List of owned components
94 RooListProxy _lowSet ; ///< Low-side variation
95 RooListProxy _highSet ; ///< High-side variation
96 RooListProxy _paramSet ; ///< interpolation parameters
97 RooListProxy _normSet ; ///< interpolation parameters
98 bool _positiveDefinite = false; ///< protect against negative and 0 bins.
99
100 std::vector<int> _interpCode;
101
102 double evaluate() const override;
103 void doEval(RooFit::EvalContext &) const override;
104
105 ClassDefOverride(PiecewiseInterpolation,4) // Sum of RooAbsReal objects
106};
107
108#endif
RooCollectionProxy< RooArgList > RooListProxy
Definition RooAbsArg.h:53
RooTemplateProxy< RooAbsReal > RooRealProxy
Compatibility typedef replacing the old RooRealProxy class.
int Int_t
Definition RtypesCore.h:45
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
char name[80]
Definition TGX11.cxx:110
RooArgList containedArgs(Action) override
const RooArgList & highList() const
bool _positiveDefinite
protect against negative and 0 bins.
const RooAbsReal * nominalHist() const
Return pointer to the nominal hist function.
RooListProxy _lowSet
Low-side variation.
RooListProxy _highSet
High-side variation.
bool isBinnedDistribution(const RooArgSet &obs) const override
WVE note: assumes nominal and alternates have identical structure, must add explicit check.
void setInterpCode(RooAbsReal &param, int code, bool silent=false)
const RooArgList & lowList() const
~PiecewiseInterpolation() override
Destructor.
TObject * clone(const char *newname) const override
RooListProxy _normSet
interpolation parameters
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
void setPositiveDefinite(bool flag=true)
RooObjCacheManager _normIntMgr
! The integration cache manager
const RooArgList & paramList() const
bool setBinIntegrator(RooArgSet &allVars)
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Advertise that all integrals can be handled internally.
RooListProxy _paramSet
interpolation parameters
const std::vector< int > & interpolationCodes() const
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
WVE note: assumes nominal and alternates have identical structure, must add explicit check.
RooArgList _ownedList
List of owned components.
RooRealProxy _nominal
The nominal value.
double evaluate() const override
Calculate and return current value of self.
void doEval(RooFit::EvalContext &) const override
Interpolate between input distributions for all values of the observable in evalData.
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Implement analytical integrations by doing appropriate weighting from component integrals functions t...
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooAbsReal()
coverity[UNINIT_CTOR] Default constructor
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
A class to maintain the context for squashing of RooFit models into code.
Implementation of a RooCacheManager<RooAbsCacheElement> that specializes in the storage of cache elem...
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Mother of all ROOT objects.
Definition TObject.h:41