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

#include <MELANCSpline_2D_fast.h>

Inheritance diagram for MELANCSpline_2D_fast:
Collaboration diagram for MELANCSpline_2D_fast:

Public Member Functions

 MELANCSpline_2D_fast ()
 
 MELANCSpline_2D_fast (const char *name, const char *title)
 
 MELANCSpline_2D_fast (const char *name, const char *title, RooAbsReal &inXVar, RooAbsReal &inYVar, const std::vector< T > &inXList, const std::vector< T > &inYList, const std::vector< std::vector< T >> &inFcnList, MELANCSplineCore::BoundaryCondition const bcBeginX_=MELANCSplineCore::bcNaturalSpline, MELANCSplineCore::BoundaryCondition const bcEndX_=MELANCSplineCore::bcNaturalSpline, MELANCSplineCore::BoundaryCondition const bcBeginY_=MELANCSplineCore::bcNaturalSpline, MELANCSplineCore::BoundaryCondition const bcEndY_=MELANCSplineCore::bcNaturalSpline, Bool_t inUseFloor=true, T inFloorEval=0, T inFloorInt=0)
 
 MELANCSpline_2D_fast (const MELANCSpline_2D_fast &other, const char *name=0)
 
virtual TObject * clone (const char *newname) const
 
virtual ~MELANCSpline_2D_fast ()
 
void setRangeValidity (const T valmin, const T valmax, const Int_t whichDirection)
 
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 ()
 
unsigned int npointsY () const
 
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) const
 
void cropValueForRange (T &val, const Int_t whichDirection) const
 
virtual std::vector< std::vector< T > > getCoefficientsPerY (const std::vector< T > &kappaX, const TMatrix_t &xAinv, const Int_t &ybin, MELANCSplineCore::BoundaryCondition const &bcBegin, MELANCSplineCore::BoundaryCondition const &bcEnd, const Int_t xbin) 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

T rangeYmin
 
T rangeYmax
 
const BoundaryCondition bcBeginX
 
const BoundaryCondition bcEndX
 
const BoundaryCondition bcBeginY
 
const BoundaryCondition bcEndY
 
RooRealProxy theYVar
 
std::vector< TYList
 
std::vector< std::vector< T > > FcnList
 
std::vector< TkappaX
 
std::vector< TkappaY
 
std::vector< std::vector< 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 11 of file MELANCSpline_2D_fast.h.

Constructor & Destructor Documentation

◆ MELANCSpline_2D_fast() [1/4]

MELANCSpline_2D_fast::MELANCSpline_2D_fast ( )

◆ MELANCSpline_2D_fast() [2/4]

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

◆ MELANCSpline_2D_fast() [3/4]

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

Definition at line 34 of file MELANCSpline_2D_fast.cc.

49  :
50  MELANCSplineCore(name, title, inXVar, inXList, inUseFloor, inFloorEval, inFloorInt),
51  rangeYmin(1), rangeYmax(-1),
52  bcBeginX(bcBeginX_), bcEndX(bcEndX_),
53  bcBeginY(bcBeginY_), bcEndY(bcEndY_),
54  theYVar("theYVar", "theYVar", this, inYVar),
55  YList(inYList),
56  FcnList(inFcnList)
57 {
58  if (npointsX()>1 && npointsY()>1){
59  // Prepare A and kappa arrays for x and y coordinates
60  int npoints;
61  Double_t det;
62 
63  vector<vector<MELANCSplineCore::T>> xA; getKappas(kappaX, 0); getAArray(kappaX, xA, bcBeginX, bcEndX);
64  npoints=kappaX.size();
65  TMatrix_t xAtrans(npoints, npoints);
66  for (int i=0; i<npoints; i++){ for (int j=0; j<npoints; j++){ xAtrans[i][j]=xA.at(i).at(j); } }
67  det=0;
68  TMatrix_t xAinv = xAtrans.Invert(&det);
69  if (det==0.){
70  coutE(InputArguments) << "MELANCSpline_2D_fast::interpolateFcn: Matrix xA could not be inverted. Something is wrong with the x coordinates of points!" << endl;
71  assert(0);
72  }
73 
74  vector<vector<MELANCSplineCore::T>> yA; getKappas(kappaY, 1); getAArray(kappaY, yA, bcBeginY, bcEndY);
75  npoints=kappaY.size();
76  TMatrix_t yAtrans(npoints, npoints);
77  for (int i=0; i<npoints; i++){ for (int j=0; j<npoints; j++){ yAtrans[i][j]=yA.at(i).at(j); } }
78  det=0;
79  TMatrix_t yAinv = yAtrans.Invert(&det);
80  if (det==0.){
81  coutE(InputArguments) << "MELANCSpline_2D_fast::interpolateFcn: Matrix yA could not be inverted. Something is wrong with the y coordinates of points!" << endl;
82  assert(0);
83  }
84 
85  // Get the grid of coefficients
86  vector<vector<vector<MELANCSplineCore::T>>> coefsAlongY; // [Ax(y),Bx(y),Cx(y),Dx(y)][xbin][ybin]
87  int npoldim=0;
88  int nxbins=0;
89  for (unsigned int j=0; j<npointsY(); j++){
90  vector<vector<MELANCSplineCore::T>> xcoefsAtYj = getCoefficientsPerY(kappaX, xAinv, j, bcBeginX, bcEndX, -1); // [ix][Ax,Bx,Cx,Dx] at each y_j
91  if (j==0){
92  nxbins=xcoefsAtYj.size();
93  npoldim=xcoefsAtYj.at(0).size();
94  for (int ipow=0; ipow<npoldim; ipow++){
95  vector<vector<MELANCSplineCore::T>> dum_xycoefarray;
96  for (int ix=0; ix<nxbins; ix++){
97  vector<MELANCSplineCore::T> dum_ycoefarray;
98  dum_xycoefarray.push_back(dum_ycoefarray);
99  }
100  coefsAlongY.push_back(dum_xycoefarray);
101  }
102  }
103  if (nxbins!=(int)xcoefsAtYj.size() || npoldim!=(int)xcoefsAtYj.at(0).size()){
104  coutE(InputArguments) << "MELANCSpline_2D_fast::interpolateFcn: nxbins!=(int)xcoefsAtYj.size() || npoldim!=(int)xcoefsAtYj.at(0).size()!" << endl;
105  assert(0);
106  }
107  for (int ix=0; ix<nxbins; ix++){
108  for (int ipow=0; ipow<npoldim; ipow++) coefsAlongY.at(ipow).at(ix).push_back(xcoefsAtYj.at(ix).at(ipow));
109  }
110  }
111 
112  for (int ix=0; ix<nxbins; ix++){
113  // Get the x coefficients interpolated across y
114  vector<vector<vector<MELANCSplineCore::T>>> xCoefs;
115  for (int ic=0; ic<npoldim; ic++){
116  vector<vector<MELANCSplineCore::T>> yCoefs = getCoefficientsAlongDirection(kappaY, yAinv, coefsAlongY.at(ic).at(ix), bcBeginY, bcEndY, -1); // [iy][A,B,C,D]
117  xCoefs.push_back(yCoefs);
118  }
119  coefficients.push_back(xCoefs);
120  }
121  }
122  else assert(0);
123 
124  RooArgSet leafset;
125  getLeafDependents(theXVar, leafset);
126  getLeafDependents(theYVar, leafset);
127  addLeafDependents(leafset);
128 
129  emptyFcnList();
130 }

◆ MELANCSpline_2D_fast() [4/4]

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

Definition at line 132 of file MELANCSpline_2D_fast.cc.

135  :
136  MELANCSplineCore(other, name),
137  rangeYmin(other.rangeYmin), rangeYmax(other.rangeYmax),
138  bcBeginX(other.bcBeginX), bcEndX(other.bcEndX),
139  bcBeginY(other.bcBeginY), bcEndY(other.bcEndY),
140  theYVar("theYVar", this, other.theYVar),
141  YList(other.YList),
142  FcnList(other.FcnList),
143  kappaX(other.kappaX),
144  kappaY(other.kappaY),
146 {}

◆ ~MELANCSpline_2D_fast()

virtual MELANCSpline_2D_fast::~MELANCSpline_2D_fast ( )
inlinevirtual

Definition at line 54 of file MELANCSpline_2D_fast.h.

54 {}

Member Function Documentation

◆ analyticalIntegral()

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

Implements MELANCSplineCore.

Definition at line 347 of file MELANCSpline_2D_fast.cc.

347  {
348  Double_t value = interpolateFcn(code, rangeName);
349  if (useFloor && value<floorInt){
350  if (verbosity>=MELANCSplineCore::kError) coutE(Integration) << "MELANCSpline_2D_fast ERROR::MELANCSpline_2D_fast(" << GetName() << ") integration returned " << value << " for code = " << code << endl;
351  value = floorInt;
352  }
353  if (verbosity==MELANCSplineCore::kVerbose){ cout << "MELANCSpline_2D_fast(" << GetName() << ")::analyticalIntegral = " << value << " for code = " << code << endl; }
354  return value;
355 }

◆ clone()

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

Implements MELANCSplineCore.

Definition at line 53 of file MELANCSpline_2D_fast.h.

53 { return new MELANCSpline_2D_fast(*this, newname); }

◆ cropValueForRange()

void MELANCSpline_2D_fast::cropValueForRange ( T val,
const Int_t  whichDirection 
) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 382 of file MELANCSpline_2D_fast.cc.

382  {
383  if (testRangeValidity(val, whichDirection)) return;
384  const T* range[2];
385  if (whichDirection==0){
386  range[0] = &rangeXmin;
387  range[1] = &rangeXmax;
388  }
389  else{
390  range[0] = &rangeYmin;
391  range[1] = &rangeYmax;
392  }
393  if (val<*(range[0])) val = *(range[0]);
394  if (val>*(range[1])) val = *(range[1]);
395 }

◆ emptyFcnList()

virtual void MELANCSpline_2D_fast::emptyFcnList ( )
inlineprotectedvirtual

Implements MELANCSplineCore.

Definition at line 62 of file MELANCSpline_2D_fast.h.

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

◆ evaluate()

Double_t MELANCSpline_2D_fast::evaluate ( ) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 313 of file MELANCSpline_2D_fast.cc.

313  {
314  Double_t value = interpolateFcn(0);
315  if (useFloor && value<floorEval){
316  if (verbosity>=MELANCSplineCore::kError) coutE(Eval) << "MELANCSpline_2D_fast ERROR::MELANCSpline_2D_fast(" << GetName() << ") evaluation returned " << value << " at (x, y) = (" << theXVar << ", " << theYVar << ")" << endl;
317  value = floorEval;
318  }
319  if (verbosity==MELANCSplineCore::kVerbose){ cout << "MELANCSpline_2D_fast(" << GetName() << ")::evaluate = " << value << " at (x, y) = (" << theXVar << ", " << theYVar << ")" << endl; }
320  return value;
321 }

◆ getAnalyticalIntegral()

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

Implements MELANCSplineCore.

Definition at line 322 of file MELANCSpline_2D_fast.cc.

322  {
323  if (_forceNumInt) return 0;
324 
325  Int_t code=1;
326 
327  RooArgSet Xdeps, Ydeps;
328  RooRealVar* rrv_x = dynamic_cast<RooRealVar*>(theXVar.absArg());
329  RooRealVar* rrv_y = dynamic_cast<RooRealVar*>(theYVar.absArg());
330  if (rrv_x==0) theXVar.absArg()->leafNodeServerList(&Xdeps, 0, true);
331  if (rrv_y==0) theYVar.absArg()->leafNodeServerList(&Ydeps, 0, true);
332 
333  if (rrv_x!=0){
334  if (Ydeps.find(*rrv_x)==0 || rrv_y!=0){
335  if (matchArgs(allVars, analVars, theXVar)) code*=2;
336  }
337  }
338  if (rrv_y!=0){
339  if (Xdeps.find(*rrv_y)==0 || rrv_x!=0){
340  if (matchArgs(allVars, analVars, theYVar)) code*=3;
341  }
342  }
343 
344  if (code==1) code=0;
345  return code;
346 }

◆ getCoefficientsPerY()

vector< vector< MELANCSplineCore::T > > MELANCSpline_2D_fast::getCoefficientsPerY ( const std::vector< T > &  kappaX,
const TMatrix_t xAinv,
const Int_t &  ybin,
MELANCSplineCore::BoundaryCondition const &  bcBegin,
MELANCSplineCore::BoundaryCondition const &  bcEnd,
const Int_t  xbin 
) const
protectedvirtual

Definition at line 306 of file MELANCSpline_2D_fast.cc.

306  {
307  vector<MELANCSplineCore::T> fcnList;
308  for (unsigned int bin=0; bin<npointsX(); bin++){ fcnList.push_back(FcnList.at(ybin).at(bin)); }
309  vector<vector<MELANCSplineCore::T>> coefs = getCoefficientsAlongDirection(kappaX, xAinv, fcnList, bcBegin, bcEnd, xbin);
310  return coefs;
311 }

◆ getKappas()

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

Implements MELANCSplineCore.

Definition at line 245 of file MELANCSpline_2D_fast.cc.

245  {
246  kappas.clear();
248 
249  Int_t npoints;
250  vector<MELANCSplineCore::T> const* coord;
251  if (whichDirection==0){
252  npoints=npointsX();
253  coord=&XList;
254  }
255  else{
256  npoints=npointsY();
257  coord=&YList;
258  }
259 
260  for (Int_t j=0; j<npoints-1; j++){
261  MELANCSplineCore::T val_j = coord->at(j);
262  MELANCSplineCore::T val_jpo = coord->at(j+1);
263  MELANCSplineCore::T val_diff = (val_jpo-val_j);
264  if (fabs(val_diff)>MELANCSplineCore::T(0)) kappa = 1./val_diff;
265  else kappa = 0;
266  kappas.push_back(kappa);
267  }
268  kappas.push_back(kappa); // Push the same kappa_(N-1)=kappa_(N-2) at the end point
269 }

◆ getTVar()

MELANCSplineCore::T MELANCSpline_2D_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 298 of file MELANCSpline_2D_fast.cc.

298  {
299  const MELANCSplineCore::T& K=kappas.at(bin);
300  vector<MELANCSplineCore::T> const* coord;
301  if (whichDirection==0) coord=&XList;
302  else coord=&YList;
303  return (val-coord->at(bin))*K;
304 }

◆ getWhichBin()

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

Implements MELANCSplineCore.

Definition at line 270 of file MELANCSpline_2D_fast.cc.

270  {
271  Int_t bin=-1;
272  MELANCSplineCore::T valj, valjpo;
273  Int_t npoints;
274  vector<MELANCSplineCore::T> const* coord;
275  if (whichDirection==0){
276  coord=&XList;
277  npoints=npointsX();
278  }
279  else{
280  coord=&YList;
281  npoints=npointsY();
282  }
283 
284  if (npoints<=1) bin=0;
285  else{
286  valjpo = coord->at(0);
287  for (Int_t j=0; j<npoints-1; j++){
288  valj = coord->at(j);
289  valjpo = coord->at(j+1);
290  if (val<valjpo && val>=valj){ bin=j; break; }
291  }
292  if (bin==-1 && val>=valjpo) bin=npoints-2;
293  else if (bin==-1) bin=0;
294  }
295 
296  return bin;
297 }

◆ interpolateFcn()

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

Implements MELANCSplineCore.

Definition at line 149 of file MELANCSpline_2D_fast.cc.

149  {
151 
152  if (verbosity==MELANCSplineCore::kVerbose){ cout << "MELANCSpline_2D_fast(" << GetName() << ")::interpolateFcn begin with code: " << code << endl; }
153 
154  // Get bins
155  Int_t xbin=-1, xbinmin=-1, xbinmax=-1, ybin=-1, ybinmin=-1, ybinmax=-1;
156  MELANCSplineCore::T tx=0, txmin=0, txmax=0, ty=0, tymin=0, tymax=0;
157  if (code==0 || code%2!=0){ // Case to just compute the value at x
158  if (!testRangeValidity(theXVar, 0)) return 0;
159  xbin = getWhichBin(theXVar, 0);
160  tx = getTVar(kappaX, theXVar, xbin, 0);
161  }
162  else{ // Case to integrate along x
163  MELANCSplineCore::T coordmin = theXVar.min(rangeName); cropValueForRange(coordmin, 0);
164  MELANCSplineCore::T coordmax = theXVar.max(rangeName); cropValueForRange(coordmax, 0);
165  xbinmin = getWhichBin(coordmin, 0);
166  txmin = getTVar(kappaX, coordmin, xbinmin, 0);
167  xbinmax = getWhichBin(coordmax, 0);
168  txmax = getTVar(kappaX, coordmax, xbinmax, 0);
169  }
170  if (code==0 || code%3!=0){ // Case to just compute the value at y
171  if (!testRangeValidity(theYVar, 1)) return 0;
172  ybin = getWhichBin(theYVar, 1);
173  ty = getTVar(kappaY, theYVar, ybin, 1);
174  }
175  else{ // Case to integrate along y
176  MELANCSplineCore::T coordmin = theYVar.min(rangeName); cropValueForRange(coordmin, 1);
177  MELANCSplineCore::T coordmax = theYVar.max(rangeName); cropValueForRange(coordmax, 1);
178  ybinmin = getWhichBin(coordmin, 1);
179  tymin = getTVar(kappaY, coordmin, ybinmin, 1);
180  ybinmax = getWhichBin(coordmax, 1);
181  tymax = getTVar(kappaY, coordmax, ybinmax, 1);
182  }
183 
184  for (int ix=0; ix<(int)coefficients.size(); ix++){
185  if (
186  (xbin>=0 && ix!=xbin)
187  ||
188  (xbinmin>=0 && xbinmax>=xbinmin && !(xbinmin<=ix && ix<=xbinmax))
189  ) continue;
190 
191  MELANCSplineCore::T txlow=0, txhigh=1;
192  if (code>0 && code%2==0){
193  if (ix==xbinmin) txlow=txmin;
194  if (ix==xbinmax) txhigh=txmax;
195  }
196  else txhigh=tx;
197 
199  if (code==0 || code%2!=0) cout << "Evaluating tx=" << txhigh << " in bin " << ix << endl;
200  else cout << "Evaluating tx[" << txlow << ", " << txhigh << "] in bin " << ix << endl;
201  }
202 
203  // Get the x coefficients interpolated across y
204  vector<MELANCSplineCore::T> xCoefs;
205  for (int ic=0; ic<(int)coefficients.at(ix).size(); ic++){
206  const vector<vector<MELANCSplineCore::T>>& yCoefs = coefficients.at(ix).at(ic);
207 
208  if (verbosity==MELANCSplineCore::kVerbose) cout << "\tCoefficient " << ic << ":\n";
209 
211  for (int iy=0; iy<(int)yCoefs.size(); iy++){
212  if (
213  (ybin>=0 && iy!=ybin)
214  ||
215  (ybinmin>=0 && ybinmax>=ybinmin && !(ybinmin<=iy && iy<=ybinmax))
216  ) continue;
217 
218  MELANCSplineCore::T tylow=0, tyhigh=1;
219  if (code>0 && code%3==0){
220  if (iy==ybinmin) tylow=tymin;
221  if (iy==ybinmax) tyhigh=tymax;
222  }
223  else tyhigh=ty;
224 
226  if (code==0 || code%3!=0) cout << "\tEvaluating ty=" << tyhigh << " in bin " << iy << endl;
227  else cout << "\tEvaluating ty[" << tylow << ", " << tyhigh << "] in bin " << iy << endl;
228  }
229 
230  theCoef += evalSplineSegment(yCoefs.at(iy), kappaY.at(iy), tyhigh, tylow, (code>0 && code%3==0));
231  }
232 
233  //if (code==0) cout << "\tCoefficient is " << theCoef << endl;
234 
235  xCoefs.push_back(theCoef);
236  }
237 
238  // Evaluate value of spline at x with coefficients evaluated at y
239  res += evalSplineSegment(xCoefs, kappaX.at(ix), txhigh, txlow, (code>0 && code%2==0));
240  }
241 
242  return res;
243 }

◆ npointsY()

unsigned int MELANCSpline_2D_fast::npointsY ( ) const
inlineprotected

Definition at line 64 of file MELANCSpline_2D_fast.h.

64 { return YList.size(); }

◆ setRangeValidity()

void MELANCSpline_2D_fast::setRangeValidity ( const T  valmin,
const T  valmax,
const Int_t  whichDirection 
)
virtual

Implements MELANCSplineCore.

Definition at line 369 of file MELANCSpline_2D_fast.cc.

369  {
370  T* range[2];
371  if (whichDirection==0){
372  range[0] = &rangeXmin;
373  range[1] = &rangeXmax;
374  }
375  else{
376  range[0] = &rangeYmin;
377  range[1] = &rangeYmax;
378  }
379  *(range[0])=valmin;
380  *(range[1])=valmax;
381 }

◆ testRangeValidity()

Bool_t MELANCSpline_2D_fast::testRangeValidity ( const T val,
const Int_t  whichDirection 
) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 357 of file MELANCSpline_2D_fast.cc.

357  {
358  const T* range[2];
359  if (whichDirection==0){
360  range[0] = &rangeXmin;
361  range[1] = &rangeXmax;
362  }
363  else{
364  range[0] = &rangeYmin;
365  range[1] = &rangeYmax;
366  }
367  return (*(range[0])>*(range[1]) || (val>=*(range[0]) && val<=*(range[1])));
368 }

Member Data Documentation

◆ bcBeginX

const BoundaryCondition MELANCSpline_2D_fast::bcBeginX
protected

Definition at line 16 of file MELANCSpline_2D_fast.h.

◆ bcBeginY

const BoundaryCondition MELANCSpline_2D_fast::bcBeginY
protected

Definition at line 18 of file MELANCSpline_2D_fast.h.

◆ bcEndX

const BoundaryCondition MELANCSpline_2D_fast::bcEndX
protected

Definition at line 17 of file MELANCSpline_2D_fast.h.

◆ bcEndY

const BoundaryCondition MELANCSpline_2D_fast::bcEndY
protected

Definition at line 19 of file MELANCSpline_2D_fast.h.

◆ coefficients

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

Definition at line 28 of file MELANCSpline_2D_fast.h.

◆ FcnList

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

Definition at line 24 of file MELANCSpline_2D_fast.h.

◆ kappaX

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

Definition at line 26 of file MELANCSpline_2D_fast.h.

◆ kappaY

std::vector<T> MELANCSpline_2D_fast::kappaY
protected

Definition at line 27 of file MELANCSpline_2D_fast.h.

◆ rangeYmax

T MELANCSpline_2D_fast::rangeYmax
protected

Definition at line 14 of file MELANCSpline_2D_fast.h.

◆ rangeYmin

T MELANCSpline_2D_fast::rangeYmin
protected

Definition at line 13 of file MELANCSpline_2D_fast.h.

◆ theYVar

RooRealProxy MELANCSpline_2D_fast::theYVar
protected

Definition at line 21 of file MELANCSpline_2D_fast.h.

◆ YList

std::vector<T> MELANCSpline_2D_fast::YList
protected

Definition at line 22 of file MELANCSpline_2D_fast.h.


The documentation for this class was generated from the following files:
MELANCSplineCore::kVerbose
@ kVerbose
Definition: MELANCSplineCore.h:25
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_2D_fast::kappaY
std::vector< T > kappaY
Definition: MELANCSpline_2D_fast.h:27
MELANCSpline_2D_fast::theYVar
RooRealProxy theYVar
Definition: MELANCSpline_2D_fast.h:21
MELANCSpline_2D_fast::rangeYmin
T rangeYmin
Definition: MELANCSpline_2D_fast.h:13
MELANCSpline_2D_fast::coefficients
std::vector< std::vector< std::vector< std::vector< T > > > > coefficients
Definition: MELANCSpline_2D_fast.h:28
TNumericUtil::KahanAccumulator
Definition: MELAAccumulators.h:27
MELANCSplineCore::rangeXmax
T rangeXmax
Definition: MELANCSplineCore.h:73
MELANCSpline_2D_fast::getTVar
virtual T getTVar(const std::vector< T > &kappas, const T &val, const Int_t &bin, const Int_t whichDirection) const
Definition: MELANCSpline_2D_fast.cc:298
MELANCSpline_2D_fast::getKappas
virtual void getKappas(std::vector< T > &kappas, const Int_t whichDirection)
Definition: MELANCSpline_2D_fast.cc:245
MELANCSpline_2D_fast::npointsY
unsigned int npointsY() const
Definition: MELANCSpline_2D_fast.h:64
MELANCSplineCore::getLeafDependents
void getLeafDependents(RooRealProxy &proxy, RooArgSet &set)
Definition: MELANCSplineCore.cc:275
MELANCSpline_2D_fast::kappaX
std::vector< T > kappaX
Definition: MELANCSpline_2D_fast.h:26
MELANCSpline_2D_fast::getWhichBin
virtual Int_t getWhichBin(const T &val, const Int_t whichDirection) const
Definition: MELANCSpline_2D_fast.cc:270
MELANCSpline_2D_fast::MELANCSpline_2D_fast
MELANCSpline_2D_fast()
MELANCSplineCore::kError
@ kError
Definition: MELANCSplineCore.h:24
MELANCSplineCore::rangeXmin
T rangeXmin
Definition: MELANCSplineCore.h:72
MELANCSpline_2D_fast::testRangeValidity
Bool_t testRangeValidity(const T &val, const Int_t whichDirection) const
Definition: MELANCSpline_2D_fast.cc:357
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_2D_fast::emptyFcnList
virtual void emptyFcnList()
Definition: MELANCSpline_2D_fast.h:62
MELANCSpline_2D_fast::getCoefficientsPerY
virtual std::vector< std::vector< T > > getCoefficientsPerY(const std::vector< T > &kappaX, const TMatrix_t &xAinv, const Int_t &ybin, MELANCSplineCore::BoundaryCondition const &bcBegin, MELANCSplineCore::BoundaryCondition const &bcEnd, const Int_t xbin) const
Definition: MELANCSpline_2D_fast.cc:306
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_2D_fast::YList
std::vector< T > YList
Definition: MELANCSpline_2D_fast.h:22
MELANCSplineCore::MELANCSplineCore
MELANCSplineCore()
MELANCSpline_2D_fast::rangeYmax
T rangeYmax
Definition: MELANCSpline_2D_fast.h:14
MELANCSplineCore::npointsX
unsigned int npointsX() const
Definition: MELANCSplineCore.h:84
MELANCSpline_2D_fast::bcEndX
const BoundaryCondition bcEndX
Definition: MELANCSpline_2D_fast.h:17
dd_global::cout
integer cout
Definition: DD_global.F90:21
MELANCSplineCore::useFloor
Bool_t useFloor
Definition: MELANCSplineCore.h:68
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
MELANCSpline_2D_fast::bcBeginY
const BoundaryCondition bcBeginY
Definition: MELANCSpline_2D_fast.h:18
MELANCSplineCore::addLeafDependents
void addLeafDependents(RooArgSet &set)
Definition: MELANCSplineCore.cc:280
MELANCSpline_2D_fast::cropValueForRange
void cropValueForRange(T &val, const Int_t whichDirection) const
Definition: MELANCSpline_2D_fast.cc:382
MELANCSplineCore::verbosity
VerbosityLevel verbosity
Definition: MELANCSplineCore.h:67
MELANCSpline_2D_fast::FcnList
std::vector< std::vector< T > > FcnList
Definition: MELANCSpline_2D_fast.h:24
MELANCSpline_2D_fast::bcBeginX
const BoundaryCondition bcBeginX
Definition: MELANCSpline_2D_fast.h:16
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
MELANCSpline_2D_fast::bcEndY
const BoundaryCondition bcEndY
Definition: MELANCSpline_2D_fast.h:19
MELANCSpline_2D_fast::interpolateFcn
virtual T interpolateFcn(Int_t code, const char *rangeName=0) const
Definition: MELANCSpline_2D_fast.cc:149