JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
MELANCSpline_1D_fast Class Reference

#include <MELANCSpline_1D_fast.h>

Inheritance diagram for MELANCSpline_1D_fast:
Collaboration diagram for MELANCSpline_1D_fast:

Public Member Functions

 MELANCSpline_1D_fast ()
 
 MELANCSpline_1D_fast (const char *name, const char *title)
 
 MELANCSpline_1D_fast (const char *name, const char *title, RooAbsReal &inXVar, const std::vector< T > &inXList, const std::vector< T > &inFcnList, MELANCSplineCore::BoundaryCondition const bcBeginX_=MELANCSplineCore::bcNaturalSpline, MELANCSplineCore::BoundaryCondition const bcEndX_=MELANCSplineCore::bcNaturalSpline, Bool_t inUseFloor=true, T inFloorEval=0, T inFloorInt=0)
 
 MELANCSpline_1D_fast (const MELANCSpline_1D_fast &other, const char *name=0)
 
virtual TObject * clone (const char *newname) const
 
virtual ~MELANCSpline_1D_fast ()
 
void setRangeValidity (const T valmin, const T valmax, const Int_t whichDirection=0)
 
virtual Int_t getAnalyticalIntegral (RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
 
virtual Double_t analyticalIntegral (Int_t code, const char *rangeName=0) const
 
- Public Member Functions inherited from MELANCSplineCore
 MELANCSplineCore ()
 
 MELANCSplineCore (const char *name, const char *title)
 
 MELANCSplineCore (const char *name, const char *title, RooAbsReal &inXVar, const std::vector< T > &inXList, Bool_t inUseFloor=true, T inFloorEval=1e-15, T inFloorInt=1e-10)
 
 MELANCSplineCore (const MELANCSplineCore &other, const char *name=0)
 
virtual ~MELANCSplineCore ()
 
virtual void setVerbosity (VerbosityLevel flag)
 
void setEvalFloor (T val)
 
void setIntFloor (T val)
 
void doFloor (Bool_t flag)
 

Protected Member Functions

virtual void emptyFcnList ()
 
virtual Int_t getWhichBin (const T &val, const Int_t whichDirection) const
 
virtual T getTVar (const std::vector< T > &kappas, const T &val, const Int_t &bin, const Int_t whichDirection) const
 
virtual void getKappas (std::vector< T > &kappas, const Int_t whichDirection)
 
Bool_t testRangeValidity (const T &val, const Int_t whichDirection=0) const
 
void cropValueForRange (T &val, const Int_t whichDirection=0) const
 
virtual T interpolateFcn (Int_t code, const char *rangeName=0) const
 
virtual Double_t evaluate () const
 
- Protected Member Functions inherited from MELANCSplineCore
void getLeafDependents (RooRealProxy &proxy, RooArgSet &set)
 
void addLeafDependents (RooArgSet &set)
 
unsigned int npointsX () const
 
virtual void getBArray (const std::vector< T > &kappas, const std::vector< T > &fcnList, std::vector< T > &BArray, BoundaryCondition const &bcBegin, BoundaryCondition const &bcEnd) const
 
virtual void getAArray (const std::vector< T > &kappas, std::vector< std::vector< T >> &AArray, BoundaryCondition const &bcBegin, BoundaryCondition const &bcEnd) const
 
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
 
virtual std::vector< TgetCoefficients (const TVector_t &S, const std::vector< T > &kappas, const std::vector< T > &fcnList, const Int_t &bin) const
 
virtual T evalSplineSegment (const std::vector< T > &coefs, const T &kappa, const T &tup, const T &tdn, Bool_t doIntegrate=false) const
 

Protected Attributes

const BoundaryCondition bcBeginX
 
const BoundaryCondition bcEndX
 
std::vector< TFcnList
 
std::vector< TkappaX
 
std::vector< std::vector< T > > coefficients
 
- Protected Attributes inherited from MELANCSplineCore
VerbosityLevel verbosity
 
Bool_t useFloor
 
T floorEval
 
T floorInt
 
T rangeXmin
 
T rangeXmax
 
RooRealProxy theXVar
 
RooListProxy leafDepsList
 
std::vector< TXList
 

Additional Inherited Members

- Public Types inherited from MELANCSplineCore
enum  VerbosityLevel { kSilent, kError, kVerbose }
 
enum  BoundaryCondition {
  bcApproximatedSlope, bcClamped, bcApproximatedSecondDerivative, bcNaturalSpline,
  bcQuadratic, bcQuadraticWithNullSlope, NBoundaryConditions
}
 
typedef Float_t T
 
typedef TMatrixT< TTMatrix_t
 
typedef TVectorT< TTVector_t
 

Detailed Description

Definition at line 12 of file MELANCSpline_1D_fast.h.

Constructor & Destructor Documentation

◆ MELANCSpline_1D_fast() [1/4]

MELANCSpline_1D_fast::MELANCSpline_1D_fast ( )

◆ MELANCSpline_1D_fast() [2/4]

MELANCSpline_1D_fast::MELANCSpline_1D_fast ( const char *  name,
const char *  title 
)

◆ MELANCSpline_1D_fast() [3/4]

MELANCSpline_1D_fast::MELANCSpline_1D_fast ( const char *  name,
const char *  title,
RooAbsReal &  inXVar,
const std::vector< T > &  inXList,
const std::vector< T > &  inFcnList,
MELANCSplineCore::BoundaryCondition const  bcBeginX_ = MELANCSplineCore::bcNaturalSpline,
MELANCSplineCore::BoundaryCondition const  bcEndX_ = MELANCSplineCore::bcNaturalSpline,
Bool_t  inUseFloor = true,
T  inFloorEval = 0,
T  inFloorInt = 0 
)

Definition at line 28 of file MELANCSpline_1D_fast.cc.

39  :
40  MELANCSplineCore(name, title, inXVar, inXList, inUseFloor, inFloorEval, inFloorInt),
41  bcBeginX(bcBeginX_), bcEndX(bcEndX_),
42  FcnList(inFcnList)
43 {
44  if (npointsX()>1){
45  int npoints;
46 
47  vector<vector<MELANCSplineCore::T>> xA; getKappas(kappaX, 0); getAArray(kappaX, xA, bcBeginX, bcEndX);
48  npoints=kappaX.size();
49  TMatrix_t xAtrans(npoints, npoints);
50  for (int i=0; i<npoints; i++){ for (int j=0; j<npoints; j++){ xAtrans[i][j]=xA.at(i).at(j); } }
51  Double_t det=0;
52  TMatrix_t xAinv = xAtrans.Invert(&det);
53  if (det==0.){
54  coutE(InputArguments) << "MELANCSpline_1D::interpolateFcn: Matrix xA could not be inverted. Something is wrong with the x coordinates of points!" << endl;
55  assert(0);
56  }
57 
59  }
60  else assert(0);
61 
62  RooArgSet leafset;
63  getLeafDependents(theXVar, leafset);
64  addLeafDependents(leafset);
65 
66  emptyFcnList();
67 }

◆ MELANCSpline_1D_fast() [4/4]

MELANCSpline_1D_fast::MELANCSpline_1D_fast ( const MELANCSpline_1D_fast other,
const char *  name = 0 
)

Definition at line 69 of file MELANCSpline_1D_fast.cc.

72  :
73  MELANCSplineCore(other, name),
74  bcBeginX(other.bcBeginX), bcEndX(other.bcEndX),
75  FcnList(other.FcnList),
76  kappaX(other.kappaX),
78 {}

◆ ~MELANCSpline_1D_fast()

virtual MELANCSpline_1D_fast::~MELANCSpline_1D_fast ( )
inlinevirtual

Definition at line 42 of file MELANCSpline_1D_fast.h.

42 {}

Member Function Documentation

◆ analyticalIntegral()

Double_t MELANCSpline_1D_fast::analyticalIntegral ( Int_t  code,
const char *  rangeName = 0 
) const
virtual

Implements MELANCSplineCore.

Definition at line 197 of file MELANCSpline_1D_fast.cc.

197  {
198  Double_t value = interpolateFcn(code, rangeName);
199  if (useFloor && value<floorInt){
200  if (verbosity>=MELANCSplineCore::kError) coutE(Integration) << "MELANCSpline_1D_fast ERROR::MELANCSpline_1D_fast(" << GetName() << ") integration returned " << value << " for code = " << code << endl;
201  value = floorInt;
202  }
203  if (verbosity==MELANCSplineCore::kVerbose){ cout << "MELANCSpline_1D_fast(" << GetName() << ")::analyticalIntegral = " << value << " for code = " << code << endl; }
204  return value;
205 }

◆ clone()

virtual TObject* MELANCSpline_1D_fast::clone ( const char *  newname) const
inlinevirtual

Implements MELANCSplineCore.

Definition at line 41 of file MELANCSpline_1D_fast.h.

41 { return new MELANCSpline_1D_fast(*this, newname); }

◆ cropValueForRange()

void MELANCSpline_1D_fast::cropValueForRange ( T val,
const Int_t  whichDirection = 0 
) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 220 of file MELANCSpline_1D_fast.cc.

220  {
221  if (testRangeValidity(val)) return;
222  const T* range[2];
223  range[0] = &rangeXmin;
224  range[1] = &rangeXmax;
225  if (val<*(range[0])) val = *(range[0]);
226  if (val>*(range[1])) val = *(range[1]);
227 }

◆ emptyFcnList()

virtual void MELANCSpline_1D_fast::emptyFcnList ( )
inlineprotectedvirtual

Implements MELANCSplineCore.

Definition at line 50 of file MELANCSpline_1D_fast.h.

50 { std::vector<T> tmp; FcnList.swap(tmp); }

◆ evaluate()

Double_t MELANCSpline_1D_fast::evaluate ( ) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 170 of file MELANCSpline_1D_fast.cc.

170  {
171  Double_t value = interpolateFcn(0);
172  if (useFloor && value<floorEval){
173  if (verbosity>=MELANCSplineCore::kError) coutE(Eval) << "MELANCSpline_1D_fast ERROR::MELANCSpline_1D_fast(" << GetName() << ") evaluation returned " << value << " at x = " << theXVar << endl;
174  value = floorEval;
175  }
177  cout << "MELANCSpline_1D_fast(" << GetName() << ")::evaluate = " << value << " at x = " << theXVar << endl;
178  RooArgSet Xdeps; theXVar.absArg()->leafNodeServerList(&Xdeps, 0, true);
179  TIterator* iter = Xdeps.createIterator();
180  RooAbsArg* var;
181  while ((var = (RooAbsArg*)iter->Next())){
182  cout << var->GetName() << " value = " << dynamic_cast<RooAbsReal*>(var)->getVal() << endl;
183  }
184  delete iter;
185  cout << endl;
186  }
187  return value;
188 }

◆ getAnalyticalIntegral()

Int_t MELANCSpline_1D_fast::getAnalyticalIntegral ( RooArgSet &  allVars,
RooArgSet &  analVars,
const char *  rangeName = 0 
) const
virtual

Implements MELANCSplineCore.

Definition at line 189 of file MELANCSpline_1D_fast.cc.

189  {
190  if (_forceNumInt) return 0;
191  Int_t code=0;
192  if (dynamic_cast<RooRealVar*>(theXVar.absArg())!=0){
193  if (matchArgs(allVars, analVars, theXVar)) code=2;
194  }
195  return code;
196 }

◆ getKappas()

void MELANCSpline_1D_fast::getKappas ( std::vector< T > &  kappas,
const Int_t  whichDirection 
)
protectedvirtual

Implements MELANCSplineCore.

Definition at line 124 of file MELANCSpline_1D_fast.cc.

124  {
125  kappas.clear();
127 
128  Int_t npoints;
129  vector<MELANCSplineCore::T> const* coord;
130  npoints=npointsX();
131  coord=&XList;
132 
133  for (Int_t j=0; j<npoints-1; j++){
134  MELANCSplineCore::T val_j = coord->at(j);
135  MELANCSplineCore::T val_jpo = coord->at(j+1);
136  MELANCSplineCore::T val_diff = (val_jpo-val_j);
137  if (fabs(val_diff)>MELANCSplineCore::T(0)) kappa = 1./val_diff;
138  else kappa = 0;
139  kappas.push_back(kappa);
140  }
141  kappas.push_back(kappa); // Push the same kappa_(N-1)=kappa_(N-2) at the end point
142 }

◆ getTVar()

MELANCSplineCore::T MELANCSpline_1D_fast::getTVar ( const std::vector< T > &  kappas,
const T val,
const Int_t &  bin,
const Int_t  whichDirection 
) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 165 of file MELANCSpline_1D_fast.cc.

165  {
166  const MELANCSplineCore::T& K=kappas.at(bin);
167  return (val-XList.at(bin))*K;
168 }

◆ getWhichBin()

Int_t MELANCSpline_1D_fast::getWhichBin ( const T val,
const Int_t  whichDirection 
) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 143 of file MELANCSpline_1D_fast.cc.

143  {
144  Int_t bin=-1;
145  MELANCSplineCore::T valj, valjpo;
146  Int_t npoints;
147  vector<MELANCSplineCore::T> const* coord;
148  coord=&XList;
149  npoints=npointsX();
150 
151  if (npoints<=1) bin=0;
152  else{
153  valjpo = coord->at(0);
154  for (Int_t j=0; j<npoints-1; j++){
155  valj = coord->at(j);
156  valjpo = coord->at(j+1);
157  if (val<valjpo && val>=valj){ bin=j; break; }
158  }
159  if (bin==-1 && val>=valjpo) bin=npoints-2;
160  else if (bin==-1) bin=0;
161  }
162 
163  return bin;
164 }

◆ interpolateFcn()

MELANCSplineCore::T MELANCSpline_1D_fast::interpolateFcn ( Int_t  code,
const char *  rangeName = 0 
) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 80 of file MELANCSpline_1D_fast.cc.

80  {
82 
83  if (verbosity==MELANCSplineCore::kVerbose) cout << "MELANCSpline_1D_fast(" << GetName() << ")::interpolateFcn begin with code: " << code << endl;
84 
85  // Get bins
86  Int_t xbin=-1, xbinmin=-1, xbinmax=-1;
87  MELANCSplineCore::T tx=0, txmin=0, txmax=0;
88  if (code==0 || code%2!=0){ // Case to just compute the value at x
89  if (!testRangeValidity(theXVar)) return 0;
90  xbin = getWhichBin(theXVar, 0);
91  tx = getTVar(kappaX, theXVar, xbin, 0);
92  }
93  else{ // Case to integrate along x
94  MELANCSplineCore::T coordmin = theXVar.min(rangeName); cropValueForRange(coordmin);
95  MELANCSplineCore::T coordmax = theXVar.max(rangeName); cropValueForRange(coordmax);
96  xbinmin = getWhichBin(coordmin, 0);
97  txmin = getTVar(kappaX, coordmin, xbinmin, 0);
98  xbinmax = getWhichBin(coordmax, 0);
99  txmax = getTVar(kappaX, coordmax, xbinmax, 0);
100  }
101 
102  int nxbins = (int)coefficients.size();
103  for (int ix=0; ix<nxbins; ix++){
104  if (
105  (xbin>=0 && ix!=xbin)
106  ||
107  (xbinmin>=0 && xbinmax>=xbinmin && !(xbinmin<=ix && ix<=xbinmax))
108  ) continue;
109 
110  MELANCSplineCore::T txlow=0, txhigh=1;
111  if (code>0 && code%2==0){
112  if (ix==xbinmin) txlow=txmin;
113  if (ix==xbinmax) txhigh=txmax;
114  }
115  else txhigh=tx;
116 
117  // Get the x coefficients at bin ix and evaluate value of spline at x
118  res += evalSplineSegment(coefficients.at(ix), kappaX.at(ix), txhigh, txlow, (code>0 && code%2==0));
119  }
120 
121  return res;
122 }

◆ setRangeValidity()

void MELANCSpline_1D_fast::setRangeValidity ( const T  valmin,
const T  valmax,
const Int_t  whichDirection = 0 
)
virtual

Implements MELANCSplineCore.

Definition at line 213 of file MELANCSpline_1D_fast.cc.

213  {
214  T* range[2];
215  range[0] = &rangeXmin;
216  range[1] = &rangeXmax;
217  *(range[0])=valmin;
218  *(range[1])=valmax;
219 }

◆ testRangeValidity()

Bool_t MELANCSpline_1D_fast::testRangeValidity ( const T val,
const Int_t  whichDirection = 0 
) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 207 of file MELANCSpline_1D_fast.cc.

207  {
208  const T* range[2];
209  range[0] = &rangeXmin;
210  range[1] = &rangeXmax;
211  return (*(range[0])>*(range[1]) || (val>=*(range[0]) && val<=*(range[1])));
212 }

Member Data Documentation

◆ bcBeginX

const BoundaryCondition MELANCSpline_1D_fast::bcBeginX
protected

Definition at line 14 of file MELANCSpline_1D_fast.h.

◆ bcEndX

const BoundaryCondition MELANCSpline_1D_fast::bcEndX
protected

Definition at line 15 of file MELANCSpline_1D_fast.h.

◆ coefficients

std::vector<std::vector<T> > MELANCSpline_1D_fast::coefficients
protected

Definition at line 20 of file MELANCSpline_1D_fast.h.

◆ FcnList

std::vector<T> MELANCSpline_1D_fast::FcnList
protected

Definition at line 17 of file MELANCSpline_1D_fast.h.

◆ kappaX

std::vector<T> MELANCSpline_1D_fast::kappaX
protected

Definition at line 19 of file MELANCSpline_1D_fast.h.


The documentation for this class was generated from the following files:
MELANCSplineCore::kVerbose
@ kVerbose
Definition: MELANCSplineCore.h:25
MELANCSpline_1D_fast::interpolateFcn
virtual T interpolateFcn(Int_t code, const char *rangeName=0) const
Definition: MELANCSpline_1D_fast.cc:80
value
pymela::gHIGGS_KAPPA value("gHIGGS_KAPPA_TILDE", pymela::gHIGGS_KAPPA_TILDE) .value("SIZE_HQQ"
MELANCSplineCore::T
Float_t T
Definition: MELANCSplineCore.h:18
MELANCSpline_1D_fast::getTVar
virtual T getTVar(const std::vector< T > &kappas, const T &val, const Int_t &bin, const Int_t whichDirection) const
Definition: MELANCSpline_1D_fast.cc:165
TNumericUtil::KahanAccumulator
Definition: MELAAccumulators.h:27
MELANCSplineCore::rangeXmax
T rangeXmax
Definition: MELANCSplineCore.h:73
MELANCSpline_1D_fast::MELANCSpline_1D_fast
MELANCSpline_1D_fast()
MELANCSplineCore::getLeafDependents
void getLeafDependents(RooRealProxy &proxy, RooArgSet &set)
Definition: MELANCSplineCore.cc:275
MELANCSpline_1D_fast::FcnList
std::vector< T > FcnList
Definition: MELANCSpline_1D_fast.h:17
MELANCSpline_1D_fast::emptyFcnList
virtual void emptyFcnList()
Definition: MELANCSpline_1D_fast.h:50
MELANCSplineCore::kError
@ kError
Definition: MELANCSplineCore.h:24
MELANCSpline_1D_fast::bcEndX
const BoundaryCondition bcEndX
Definition: MELANCSpline_1D_fast.h:15
MELANCSplineCore::rangeXmin
T rangeXmin
Definition: MELANCSplineCore.h:72
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
MELANCSpline_1D_fast::getWhichBin
virtual Int_t getWhichBin(const T &val, const Int_t whichDirection) const
Definition: MELANCSpline_1D_fast.cc:143
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::TMatrix_t
TMatrixT< T > TMatrix_t
Definition: MELANCSplineCore.h:19
MELANCSpline_1D_fast::getKappas
virtual void getKappas(std::vector< T > &kappas, const Int_t whichDirection)
Definition: MELANCSpline_1D_fast.cc:124
MELANCSplineCore::MELANCSplineCore
MELANCSplineCore()
MELANCSplineCore::npointsX
unsigned int npointsX() const
Definition: MELANCSplineCore.h:84
MELANCSpline_1D_fast::cropValueForRange
void cropValueForRange(T &val, const Int_t whichDirection=0) const
Definition: MELANCSpline_1D_fast.cc:220
MELANCSpline_1D_fast::testRangeValidity
Bool_t testRangeValidity(const T &val, const Int_t whichDirection=0) const
Definition: MELANCSpline_1D_fast.cc:207
MELANCSpline_1D_fast::bcBeginX
const BoundaryCondition bcBeginX
Definition: MELANCSpline_1D_fast.h:14
dd_global::cout
integer cout
Definition: DD_global.F90:21
MELANCSplineCore::useFloor
Bool_t useFloor
Definition: MELANCSplineCore.h:68
MELANCSpline_1D_fast::kappaX
std::vector< T > kappaX
Definition: MELANCSpline_1D_fast.h:19
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::addLeafDependents
void addLeafDependents(RooArgSet &set)
Definition: MELANCSplineCore.cc:280
MELANCSpline_1D_fast::coefficients
std::vector< std::vector< T > > coefficients
Definition: MELANCSpline_1D_fast.h:20
MELANCSplineCore::verbosity
VerbosityLevel verbosity
Definition: MELANCSplineCore.h:67
MELANCSplineCore::theXVar
RooRealProxy theXVar
Definition: MELANCSplineCore.h:75
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