JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
MELANCSplineCore.h
Go to the documentation of this file.
1 #ifndef MELANCSPLINECORE
2 #define MELANCSPLINECORE
3 
4 #include <vector>
5 #include "MELAAccumulators.h"
6 #include "TMatrixD.h"
7 #include "TVectorD.h"
8 #include "RooAbsReal.h"
9 #include "RooRealVar.h"
10 #include "RooRealProxy.h"
11 #include "RooListProxy.h"
12 #include "RooConstVar.h"
13 #include "RooArgList.h"
14 #include "RooMsgService.h"
15 
16 class MELANCSplineCore : public RooAbsReal{
17 public:
18  typedef Float_t T;
19  typedef TMatrixT<T> TMatrix_t;
20  typedef TVectorT<T> TVector_t;
21 
26  };
27 
29  bcApproximatedSlope, // D1 is the same as deltaY/deltaX in the first/last segment
30  bcClamped, // D1=0 at endpoint
31  bcApproximatedSecondDerivative, // D2 is approximated
32  bcNaturalSpline, // D2=0 at endpoint
33  bcQuadratic, // Coefficient D=0
34  bcQuadraticWithNullSlope, // Coefficient D=0 and D1=0 at endpoint
36  };
37 
40  const char* name,
41  const char* title
42  );
44  const char* name,
45  const char* title,
46  RooAbsReal& inXVar,
47  const std::vector<T>& inXList,
48  Bool_t inUseFloor=true,
49  T inFloorEval=1e-15,
50  T inFloorInt=1e-10
51  );
52  MELANCSplineCore(const MELANCSplineCore& other, const char* name=0);
53  virtual TObject* clone(const char* newname)const = 0;
54  inline virtual ~MELANCSplineCore(){}
55 
56  virtual void setVerbosity(VerbosityLevel flag);
57  void setEvalFloor(T val);
58  void setIntFloor(T val);
59  void doFloor(Bool_t flag);
60 
61  virtual void setRangeValidity(const T valmin, const T valmax, const Int_t whichDirection) = 0;
62 
63  virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0)const = 0;
64  virtual Double_t analyticalIntegral(Int_t code, const char* rangeName=0)const = 0;
65 
66 protected:
68  Bool_t useFloor;
71 
74 
75  RooRealProxy theXVar;
76  RooListProxy leafDepsList;
77  std::vector<T> XList;
78 
79  virtual void emptyFcnList() = 0;
80 
81  void getLeafDependents(RooRealProxy& proxy, RooArgSet& set);
82  void addLeafDependents(RooArgSet& set);
83 
84  unsigned int npointsX()const{ return XList.size(); }
85 
86  virtual Int_t getWhichBin(const T& val, const Int_t whichDirection)const = 0;
87  virtual void getKappas(std::vector<T>& kappas, const Int_t whichDirection) = 0;
88  virtual T getTVar(const std::vector<T>& kappas, const T& val, const Int_t& bin, const Int_t whichDirection)const = 0;
89 
90  virtual Bool_t testRangeValidity(const T& val, const Int_t whichDirection)const = 0;
91  virtual void cropValueForRange(T& val, const Int_t whichDirection)const = 0;
92 
93  virtual T interpolateFcn(Int_t code, const char* rangeName=0)const = 0;
94  virtual Double_t evaluate()const = 0;
95 
96  virtual void getBArray(const std::vector<T>& kappas, const std::vector<T>& fcnList, std::vector<T>& BArray, BoundaryCondition const& bcBegin, BoundaryCondition const& bcEnd)const;
97  virtual void getAArray(const std::vector<T>& kappas, std::vector<std::vector<T>>& AArray, BoundaryCondition const& bcBegin, BoundaryCondition const& bcEnd)const;
98  virtual std::vector<std::vector<T>> getCoefficientsAlongDirection(const std::vector<T>& kappas, const TMatrix_t& Ainv, const std::vector<T>& fcnList, BoundaryCondition const& bcBegin, BoundaryCondition const& bcEnd, const Int_t pickBin)const;
99  virtual std::vector<T> getCoefficients(const TVector_t& S, const std::vector<T>& kappas, const std::vector<T>& fcnList, const Int_t& bin)const;
100 
101  virtual T evalSplineSegment(const std::vector<T>& coefs, const T& kappa, const T& tup, const T& tdn, Bool_t doIntegrate=false)const;
102 
103 
104  ClassDef(MELANCSplineCore, 0)
105 
106 };
107 
108 #endif
MELANCSplineCore::testRangeValidity
virtual Bool_t testRangeValidity(const T &val, const Int_t whichDirection) const =0
MELANCSplineCore::kVerbose
@ kVerbose
Definition: MELANCSplineCore.h:25
MELANCSplineCore::T
Float_t T
Definition: MELANCSplineCore.h:18
MELANCSplineCore::VerbosityLevel
VerbosityLevel
Definition: MELANCSplineCore.h:22
MELANCSplineCore::getAnalyticalIntegral
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const =0
MELANCSplineCore::NBoundaryConditions
@ NBoundaryConditions
Definition: MELANCSplineCore.h:35
MELANCSplineCore::setIntFloor
void setIntFloor(T val)
Definition: MELANCSplineCore.cc:69
MELANCSplineCore::rangeXmax
T rangeXmax
Definition: MELANCSplineCore.h:73
MELANCSplineCore::setEvalFloor
void setEvalFloor(T val)
Definition: MELANCSplineCore.cc:68
MELANCSplineCore::getLeafDependents
void getLeafDependents(RooRealProxy &proxy, RooArgSet &set)
Definition: MELANCSplineCore.cc:275
MELANCSplineCore::getWhichBin
virtual Int_t getWhichBin(const T &val, const Int_t whichDirection) const =0
MELANCSplineCore::BoundaryCondition
BoundaryCondition
Definition: MELANCSplineCore.h:28
MELANCSplineCore::~MELANCSplineCore
virtual ~MELANCSplineCore()
Definition: MELANCSplineCore.h:54
MELANCSplineCore::clone
virtual TObject * clone(const char *newname) const =0
MELANCSplineCore::kError
@ kError
Definition: MELANCSplineCore.h:24
MELANCSplineCore::emptyFcnList
virtual void emptyFcnList()=0
MELANCSplineCore::bcQuadraticWithNullSlope
@ bcQuadraticWithNullSlope
Definition: MELANCSplineCore.h:34
MELANCSplineCore::setVerbosity
virtual void setVerbosity(VerbosityLevel flag)
Definition: MELANCSplineCore.cc:67
MELANCSplineCore::rangeXmin
T rangeXmin
Definition: MELANCSplineCore.h:72
set
set(CMAKE_MACOSX_RPATH 1) set(CMAKE_CACHEFILE_DIR CMakeFiles/) cmake_minimum_required(VERSION 2.8) project(collier) enable_language(Fortran) get_filename_component(Fortran_COMPILER_NAME $
Definition: CMakeLists.txt:65
MELANCSplineCore::getTVar
virtual T getTVar(const std::vector< T > &kappas, const T &val, const Int_t &bin, const Int_t whichDirection) const =0
modparameters::kappa
complex(8), public kappa
Definition: mod_Parameters.F90:882
MELANCSplineCore::evalSplineSegment
virtual T evalSplineSegment(const std::vector< T > &coefs, const T &kappa, const T &tup, const T &tdn, Bool_t doIntegrate=false) const
Definition: MELANCSplineCore.cc:265
MELANCSplineCore::setRangeValidity
virtual void setRangeValidity(const T valmin, const T valmax, const Int_t whichDirection)=0
MELANCSplineCore::getCoefficientsAlongDirection
virtual std::vector< std::vector< T > > getCoefficientsAlongDirection(const std::vector< T > &kappas, const TMatrix_t &Ainv, const std::vector< T > &fcnList, BoundaryCondition const &bcBegin, BoundaryCondition const &bcEnd, const Int_t pickBin) const
Definition: MELANCSplineCore.cc:247
MELANCSplineCore::bcQuadratic
@ bcQuadratic
Definition: MELANCSplineCore.h:33
MELANCSplineCore::TMatrix_t
TMatrixT< T > TMatrix_t
Definition: MELANCSplineCore.h:19
MELANCSplineCore::evaluate
virtual Double_t evaluate() const =0
MELANCSplineCore::cropValueForRange
virtual void cropValueForRange(T &val, const Int_t whichDirection) const =0
MELANCSplineCore::MELANCSplineCore
MELANCSplineCore()
MELANCSplineCore::npointsX
unsigned int npointsX() const
Definition: MELANCSplineCore.h:84
MELANCSplineCore
Definition: MELANCSplineCore.h:16
MELANCSplineCore::kSilent
@ kSilent
Definition: MELANCSplineCore.h:23
MELANCSplineCore::bcApproximatedSlope
@ bcApproximatedSlope
Definition: MELANCSplineCore.h:29
MELANCSplineCore::bcClamped
@ bcClamped
Definition: MELANCSplineCore.h:30
MELANCSplineCore::useFloor
Bool_t useFloor
Definition: MELANCSplineCore.h:68
MELANCSplineCore::interpolateFcn
virtual T interpolateFcn(Int_t code, const char *rangeName=0) const =0
MELANCSplineCore::getKappas
virtual void getKappas(std::vector< T > &kappas, const Int_t whichDirection)=0
MELANCSplineCore::bcNaturalSpline
@ bcNaturalSpline
Definition: MELANCSplineCore.h:32
MELANCSplineCore::floorInt
T floorInt
Definition: MELANCSplineCore.h:70
MELANCSplineCore::XList
std::vector< T > XList
Definition: MELANCSplineCore.h:77
MELANCSplineCore::floorEval
T floorEval
Definition: MELANCSplineCore.h:69
MELANCSplineCore::doFloor
void doFloor(Bool_t flag)
Definition: MELANCSplineCore.cc:70
MELANCSplineCore::addLeafDependents
void addLeafDependents(RooArgSet &set)
Definition: MELANCSplineCore.cc:280
MELANCSplineCore::verbosity
VerbosityLevel verbosity
Definition: MELANCSplineCore.h:67
MELANCSplineCore::TVector_t
TVectorT< T > TVector_t
Definition: MELANCSplineCore.h:20
MELANCSplineCore::theXVar
RooRealProxy theXVar
Definition: MELANCSplineCore.h:75
MELANCSplineCore::getBArray
virtual void getBArray(const std::vector< T > &kappas, const std::vector< T > &fcnList, std::vector< T > &BArray, BoundaryCondition const &bcBegin, BoundaryCondition const &bcEnd) const
Definition: MELANCSplineCore.cc:72
MELAAccumulators.h
MELANCSplineCore::analyticalIntegral
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const =0
MELANCSplineCore::getAArray
virtual void getAArray(const std::vector< T > &kappas, std::vector< std::vector< T >> &AArray, BoundaryCondition const &bcBegin, BoundaryCondition const &bcEnd) const
Definition: MELANCSplineCore.cc:159
MELANCSplineCore::bcApproximatedSecondDerivative
@ bcApproximatedSecondDerivative
Definition: MELANCSplineCore.h:31
MELANCSplineCore::getCoefficients
virtual std::vector< T > getCoefficients(const TVector_t &S, const std::vector< T > &kappas, const std::vector< T > &fcnList, const Int_t &bin) const
Definition: MELANCSplineCore.cc:219
MELANCSplineCore::leafDepsList
RooListProxy leafDepsList
Definition: MELANCSplineCore.h:76