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, const RooArgList& highSet, const RooArgList& paramSet
35#ifndef ROOFIT_MEMORY_SAFE_INTERFACES
36 , bool takeOwnership=false
37#endif
38 );
39 ~PiecewiseInterpolation() override ;
40
41 PiecewiseInterpolation(const PiecewiseInterpolation& other, const char *name = nullptr);
42 TObject* clone(const char* newname) const override { return new PiecewiseInterpolation(*this, newname); }
43
44 /// Return pointer to the nominal hist function.
45 const RooAbsReal* nominalHist() const {
46 return &_nominal.arg();
47 }
48
49 // virtual double defaultErrorLevel() const ;
50
51 // void printMetaArgs(std::ostream& os) const ;
52
53 const RooArgList& lowList() const { return _lowSet ; }
54 const RooArgList& highList() const { return _highSet ; }
55 const RooArgList& paramList() const { return _paramSet ; }
56 const std::vector<int>& interpolationCodes() const { return _interpCode; }
57
58 //virtual bool forceAnalyticalInt(const RooAbsArg&) const { return true ; }
59 bool setBinIntegrator(RooArgSet& allVars) ;
60
61 Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet,const char* rangeName=nullptr) const override ;
62 double analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=nullptr) const override ;
63
64 void setPositiveDefinite(bool flag=true){_positiveDefinite=flag;}
65 bool positiveDefinite() const {return _positiveDefinite;}
66
67 void setInterpCode(RooAbsReal& param, int code, bool silent=false);
68 void setAllInterpCodes(int code);
70
71 std::list<double>* binBoundaries(RooAbsRealLValue& /*obs*/, double /*xlo*/, double /*xhi*/) const override ;
72 std::list<double>* plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const override ;
73 bool isBinnedDistribution(const RooArgSet& obs) const override ;
74
75protected:
76
78 public:
79 CacheElem() {} ;
80 ~CacheElem() override {} ;
83 ret.add(_lowIntList);
84 ret.add(_highIntList);
85 return ret ;
86 }
90 // will want std::vector<RooRealVar*> for low and high also
91 } ;
92 mutable RooObjCacheManager _normIntMgr ; ///<! The integration cache manager
93
94 RooRealProxy _nominal; ///< The nominal value
95 RooArgList _ownedList ; ///< List of owned components
96 RooListProxy _lowSet ; ///< Low-side variation
97 RooListProxy _highSet ; ///< High-side variation
98 RooListProxy _paramSet ; ///< interpolation parameters
99 RooListProxy _normSet ; ///< interpolation parameters
100 bool _positiveDefinite; ///< protect against negative and 0 bins.
101
102 std::vector<int> _interpCode;
103
104 double evaluate() const override;
105 void computeBatch(cudaStream_t*, double* output, size_t size, RooFit::Detail::DataMap const&) const override;
106
107 ClassDefOverride(PiecewiseInterpolation,4) // Sum of RooAbsReal objects
108};
109
110#endif
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
char name[80]
Definition TGX11.cxx:110
RooArgList containedArgs(Action) override
The PiecewiseInterpolation is a class that can morph distributions into each other,...
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 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.
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...
void computeBatch(cudaStream_t *, double *output, size_t size, RooFit::Detail::DataMap const &) const override
Interpolate between input distributions for all values of the observable in evalData.
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
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
Class RooObjCacheManager is an implementation of class RooCacheManager<RooAbsCacheElement> and specia...
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
const T & arg() const
Return reference to object held in proxy.
Mother of all ROOT objects.
Definition TObject.h:41
static void output()