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, Bool_t takeOwnerShip=kFALSE) ;
35 virtual ~PiecewiseInterpolation() ;
36
37 PiecewiseInterpolation(const PiecewiseInterpolation& other, const char* name = 0);
38 virtual TObject* clone(const char* newname) const { return new PiecewiseInterpolation(*this, newname); }
39
40 /// Return pointer to the nominal hist function.
41 const RooAbsReal* nominalHist() const {
42 return &_nominal.arg();
43 }
44
45 // virtual Double_t defaultErrorLevel() const ;
46
47 // void printMetaArgs(std::ostream& os) const ;
48
49 const RooArgList& lowList() const { return _lowSet ; }
50 const RooArgList& highList() const { return _highSet ; }
51 const RooArgList& paramList() const { return _paramSet ; }
52 const std::vector<int>& interpolationCodes() const { return _interpCode; }
53
54 //virtual Bool_t forceAnalyticalInt(const RooAbsArg&) const { return kTRUE ; }
56
57 Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet,const char* rangeName=0) const ;
58 Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
59
60 void setPositiveDefinite(bool flag=true){_positiveDefinite=flag;}
61
62 void setInterpCode(RooAbsReal& param, int code, bool silent=false);
63 void setAllInterpCodes(int code);
65
66 virtual std::list<Double_t>* binBoundaries(RooAbsRealLValue& /*obs*/, Double_t /*xlo*/, Double_t /*xhi*/) const ;
67 virtual std::list<Double_t>* plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const ;
68 virtual Bool_t isBinnedDistribution(const RooArgSet& obs) const ;
69
70protected:
71
73 public:
74 CacheElem() {} ;
75 virtual ~CacheElem() {} ;
78 ret.add(_lowIntList);
79 ret.add(_highIntList);
80 return ret ;
81 }
85 // will want std::vector<RooRealVar*> for low and high also
86 } ;
87 mutable RooObjCacheManager _normIntMgr ; ///<! The integration cache manager
88
89 RooRealProxy _nominal; ///< The nominal value
90 RooArgList _ownedList ; ///< List of owned components
91 RooListProxy _lowSet ; ///< Low-side variation
92 RooListProxy _highSet ; ///< High-side variation
93 RooListProxy _paramSet ; ///< interpolation parameters
94 RooListProxy _normSet ; ///< interpolation parameters
95 Bool_t _positiveDefinite; ///< protect against negative and 0 bins.
96
97 std::vector<int> _interpCode;
98
99 double evaluate() const;
100 void computeBatch(cudaStream_t*, double* output, size_t size, RooFit::Detail::DataMap const&) const;
101
102 ClassDef(PiecewiseInterpolation,4) // Sum of RooAbsReal objects
103};
104
105#endif
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
const Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
#define ClassDef(name, id)
Definition Rtypes.h:325
char name[80]
Definition TGX11.cxx:110
virtual RooArgList containedArgs(Action)
The PiecewiseInterpolation is a class that can morph distributions into each other,...
const RooArgList & highList() const
const RooAbsReal * nominalHist() const
Return pointer to the nominal hist function.
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
WVE note: assumes nominal and alternates have identical structure, must add explicit check.
RooListProxy _lowSet
Low-side variation.
Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Implement analytical integrations by doing appropriate weighting from component integrals functions t...
RooListProxy _highSet
High-side variation.
void setInterpCode(RooAbsReal &param, int code, bool silent=false)
const RooArgList & lowList() const
Bool_t _positiveDefinite
protect against negative and 0 bins.
Bool_t setBinIntegrator(RooArgSet &allVars)
RooListProxy _normSet
interpolation parameters
void setPositiveDefinite(bool flag=true)
RooObjCacheManager _normIntMgr
! The integration cache manager
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=0) const
Advertise that all integrals can be handled internally.
const RooArgList & paramList() const
double evaluate() const
Calculate and return current value of self.
void computeBatch(cudaStream_t *, double *output, size_t size, RooFit::Detail::DataMap const &) const
Interpolate between input distributions for all values of the observable in evalData.
RooListProxy _paramSet
interpolation parameters
const std::vector< int > & interpolationCodes() const
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
virtual TObject * clone(const char *newname) const
virtual ~PiecewiseInterpolation()
Destructor.
RooArgList _ownedList
List of owned components.
RooRealProxy _nominal
The nominal value.
virtual Bool_t isBinnedDistribution(const RooArgSet &obs) const
WVE note: assumes nominal and alternates have identical structure, must add explicit check.
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
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:64
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:35
RooListProxy is the concrete proxy for RooArgList objects.
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:39
const T & arg() const
Return reference to object held in proxy.
Mother of all ROOT objects.
Definition TObject.h:41
static void output(int code)
Definition gifencode.c:226