JHUGen MELA  v2.4.1
Matrix element calculations as used in JHUGen. MELA is an important tool that was used for the Higgs boson discovery and for precise measurements of its structure and interactions. Please see the website https://spin.pha.jhu.edu/ and papers cited there for more details, and kindly cite those papers when using this code.
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
MELANCSpline_3D_fast Class Reference

#include <MELANCSpline_3D_fast.h>

Inheritance diagram for MELANCSpline_3D_fast:
Inheritance graph
[legend]
Collaboration diagram for MELANCSpline_3D_fast:
Collaboration graph
[legend]

Public Member Functions

 MELANCSpline_3D_fast ()
 
 MELANCSpline_3D_fast (const char *name, const char *title)
 
 MELANCSpline_3D_fast (const char *name, const char *title, RooAbsReal &inXVar, RooAbsReal &inYVar, RooAbsReal &inZVar, const std::vector< T > &inXList, const std::vector< T > &inYList, const std::vector< T > &inZList, const std::vector< 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, MELANCSplineCore::BoundaryCondition const bcBeginZ_=MELANCSplineCore::bcNaturalSpline, MELANCSplineCore::BoundaryCondition const bcEndZ_=MELANCSplineCore::bcNaturalSpline, Bool_t inUseFloor=true, T inFloorEval=0, T inFloorInt=0)
 
 MELANCSpline_3D_fast (const MELANCSpline_3D_fast &other, const char *name=0)
 
virtual TObject * clone (const char *newname) const
 
virtual ~MELANCSpline_3D_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
 
unsigned int npointsZ () 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 > > getCoefficientsPerYPerZ (const std::vector< T > &kappaX, const TMatrix_t &xAinv, const Int_t &ybin, const Int_t &zbin, 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
 
T rangeZmin
 
T rangeZmax
 
const BoundaryCondition bcBeginX
 
const BoundaryCondition bcEndX
 
const BoundaryCondition bcBeginY
 
const BoundaryCondition bcEndY
 
const BoundaryCondition bcBeginZ
 
const BoundaryCondition bcEndZ
 
RooRealProxy theYVar
 
RooRealProxy theZVar
 
std::vector< TYList
 
std::vector< TZList
 
std::vector< std::vector< std::vector< T > > > FcnList
 
std::vector< TkappaX
 
std::vector< TkappaY
 
std::vector< TkappaZ
 
std::vector< std::vector< 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_3D_fast.h.

Constructor & Destructor Documentation

◆ MELANCSpline_3D_fast() [1/4]

MELANCSpline_3D_fast::MELANCSpline_3D_fast ( )

◆ MELANCSpline_3D_fast() [2/4]

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

◆ MELANCSpline_3D_fast() [3/4]

MELANCSpline_3D_fast::MELANCSpline_3D_fast ( const char *  name,
const char *  title,
RooAbsReal &  inXVar,
RooAbsReal &  inYVar,
RooAbsReal &  inZVar,
const std::vector< T > &  inXList,
const std::vector< T > &  inYList,
const std::vector< T > &  inZList,
const std::vector< 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,
MELANCSplineCore::BoundaryCondition const  bcBeginZ_ = MELANCSplineCore::bcNaturalSpline,
MELANCSplineCore::BoundaryCondition const  bcEndZ_ = MELANCSplineCore::bcNaturalSpline,
Bool_t  inUseFloor = true,
T  inFloorEval = 0,
T  inFloorInt = 0 
)

Definition at line 40 of file MELANCSpline_3D_fast.cc.

59  :
60  MELANCSplineCore(name, title, inXVar, inXList, inUseFloor, inFloorEval, inFloorInt),
61  rangeYmin(1), rangeYmax(-1),
62  rangeZmin(1), rangeZmax(-1),
63  bcBeginX(bcBeginX_), bcEndX(bcEndX_),
64  bcBeginY(bcBeginY_), bcEndY(bcEndY_),
65  bcBeginZ(bcBeginZ_), bcEndZ(bcEndZ_),
66  theYVar("theYVar", "theYVar", this, inYVar),
67  theZVar("theZVar", "theZVar", this, inZVar),
68  YList(inYList),
69  ZList(inZList),
70  FcnList(inFcnList)
71 {
72  if (npointsX()>1 && npointsY()>1 && npointsZ()>1){
73  // Prepare A and kappa arrays for x, y and z coordinates
74  int npoints;
75  Double_t det;
76 
77  vector<vector<MELANCSplineCore::T>> xA; getKappas(kappaX, 0); getAArray(kappaX, xA, bcBeginX, bcEndX);
78  npoints=kappaX.size();
79  TMatrix_t xAtrans(npoints, npoints);
80  for (int i=0; i<npoints; i++){ for (int j=0; j<npoints; j++){ xAtrans[i][j]=xA.at(i).at(j); } }
81  det=0;
82  TMatrix_t xAinv = xAtrans.Invert(&det);
83  if (det==0.){
84  coutE(InputArguments) << "MELANCSpline_3D_fast::interpolateFcn: Matrix xA could not be inverted. Something is wrong with the x coordinates of points!" << endl;
85  assert(0);
86  }
87 
88  vector<vector<MELANCSplineCore::T>> yA; getKappas(kappaY, 1); getAArray(kappaY, yA, bcBeginY, bcEndY);
89  npoints=kappaY.size();
90  TMatrix_t yAtrans(npoints, npoints);
91  for (int i=0; i<npoints; i++){ for (int j=0; j<npoints; j++){ yAtrans[i][j]=yA.at(i).at(j); } }
92  det=0;
93  TMatrix_t yAinv = yAtrans.Invert(&det);
94  if (det==0.){
95  coutE(InputArguments) << "MELANCSpline_3D_fast::interpolateFcn: Matrix yA could not be inverted. Something is wrong with the y coordinates of points!" << endl;
96  assert(0);
97  }
98 
99  vector<vector<MELANCSplineCore::T>> zA; getKappas(kappaZ, 2); getAArray(kappaZ, zA, bcBeginZ, bcEndZ);
100  npoints=kappaZ.size();
101  TMatrix_t zAtrans(npoints, npoints);
102  for (int i=0; i<npoints; i++){ for (int j=0; j<npoints; j++){ zAtrans[i][j]=zA.at(i).at(j); } }
103  det=0;
104  TMatrix_t zAinv = zAtrans.Invert(&det);
105  if (det==0.){
106  coutE(InputArguments) << "MELANCSpline_3D_fast::interpolateFcn: Matrix zA could not be inverted. Something is wrong with the z coordinates of points!" << endl;
107  assert(0);
108  }
109 
110  //cout << "MELANCSpline_3D_fast::MELANCSpline_3D_fast: Initial kappa arrays are set up" << endl;
111 
112  // Get the grid of coefficients
113  int nxpoldim=0;
114  int nxbins=0;
115  int nybins=0;
116  int nypoldim=0;
117  std::vector<
118  std::vector<std::vector<
119  std::vector<std::vector<MELANCSplineCore::T>>
120  >>
121  > coefsAlongZ; // [ix][A_x,B_x,C_x,D_x][iy][A_x_y,B_x_y,C_x_y,D_x_y][z_k] at each z_k
122  for (unsigned int k=0; k<npointsZ(); k++){
123  //cout << "Finding coefficients along z line " << k << endl;
124 
125  std::vector<std::vector<
126  std::vector<std::vector<MELANCSplineCore::T>>
127  >> coefficients_perZ; // [ix][A_x,B_x,C_x,D_x][iy][A_x_y,B_x_y,C_x_y,D_x_y] at each z_k
128 
129  vector<vector<vector<MELANCSplineCore::T>>> coefsAlongY; // [xbin][Ax(y),Bx(y),Cx(y),Dx(y)][ybin] in each z
130  for (unsigned int j=0; j<npointsY(); j++){
131  vector<vector<MELANCSplineCore::T>> xcoefsAtYjZk = getCoefficientsPerYPerZ(kappaX, xAinv, j, k, bcBeginX, bcEndX, -1); // [ix][Ax,Bx,Cx,Dx] at each y_j z_k
132  //cout << "\tCoefficients in y line " << j << " are found" << endl;
133  if (j==0){
134  if (k==0){
135  nxbins=xcoefsAtYjZk.size();
136  nxpoldim=xcoefsAtYjZk.at(0).size();
137  }
138  for (int ix=0; ix<nxbins; ix++){
139  vector<vector<MELANCSplineCore::T>> dum_xycoefarray;
140  for (int icx=0; icx<nxpoldim; icx++){
141  vector<MELANCSplineCore::T> dum_ycoefarray;
142  dum_xycoefarray.push_back(dum_ycoefarray);
143  }
144  coefsAlongY.push_back(dum_xycoefarray);
145  }
146  }
147  if (nxbins!=(int)xcoefsAtYjZk.size() || nxpoldim!=(int)xcoefsAtYjZk.at(0).size()){
148  coutE(InputArguments) << "MELANCSpline_3D_fast::interpolateFcn: nxbins!=(int)xcoefsAtYjZk.size() || nxpoldim!=(int)xcoefsAtYjZk.at(0).size()!" << endl;
149  assert(0);
150  }
151  for (int ix=0; ix<nxbins; ix++){
152  for (int icx=0; icx<nxpoldim; icx++) coefsAlongY.at(ix).at(icx).push_back(xcoefsAtYjZk.at(ix).at(icx));
153  }
154  }
155 
156  //cout << "\tCoefficients in each y line are found" << endl;
157 
158  for (int ix=0; ix<nxbins; ix++){
159  // Get the x coefficients interpolated across y
160  vector<vector<vector<MELANCSplineCore::T>>> xCoefs;
161  for (int icx=0; icx<nxpoldim; icx++){
162  vector<vector<MELANCSplineCore::T>> yCoefs = getCoefficientsAlongDirection(kappaY, yAinv, coefsAlongY.at(ix).at(icx), bcBeginY, bcEndY, -1); // [iy][A,B,C,D]
163  xCoefs.push_back(yCoefs);
164  }
165  coefficients_perZ.push_back(xCoefs);
166  }
167 
168  if (k==0){
169  nybins = coefficients_perZ.at(0).at(0).size();
170  nypoldim = coefficients_perZ.at(0).at(0).at(0).size();
171  for (int ix=0; ix<nxbins; ix++){
172  vector<vector<vector<vector<MELANCSplineCore::T>>>> xCoefs;
173  for (int icx=0; icx<nxpoldim; icx++){
174  vector<vector<vector<MELANCSplineCore::T>>> xCoefsAlongY;
175  for (int iy=0; iy<nybins; iy++){
176  vector<vector<MELANCSplineCore::T>> yCoefs;
177  for (int icy=0; icy<nypoldim; icy++){
178  vector<MELANCSplineCore::T> yCoefAtZj;
179  yCoefs.push_back(yCoefAtZj);
180  }
181  xCoefsAlongY.push_back(yCoefs);
182  }
183  xCoefs.push_back(xCoefsAlongY);
184  }
185  coefsAlongZ.push_back(xCoefs);
186  }
187  }
188  for (int ix=0; ix<nxbins; ix++){
189  for (int icx=0; icx<nxpoldim; icx++){
190  for (int iy=0; iy<nybins; iy++){
191  for (int icy=0; icy<nypoldim; icy++){
192  coefsAlongZ.at(ix).at(icx).at(iy).at(icy).push_back(coefficients_perZ.at(ix).at(icx).at(iy).at(icy));
193  }
194  }
195  }
196  }
197  }
198 
199  //cout << "Finding the final coefficients" << endl;
200  for (int ix=0; ix<nxbins; ix++){
201  vector<vector<vector<vector<vector<MELANCSplineCore::T>>>>> xCoefs;
202  for (int icx=0; icx<nxpoldim; icx++){
203  vector<vector<vector<vector<MELANCSplineCore::T>>>> xCoefsAlongY;
204  for (int iy=0; iy<nybins; iy++){
205  vector<vector<vector<MELANCSplineCore::T>>> yCoefs;
206  for (int icy=0; icy<nypoldim; icy++){
207  vector<vector<MELANCSplineCore::T>> yCoefsAlongZ = getCoefficientsAlongDirection(kappaZ, zAinv, coefsAlongZ.at(ix).at(icx).at(iy).at(icy), bcBeginZ, bcEndZ, -1); // [iz][A,B,C,D]
208  yCoefs.push_back(yCoefsAlongZ);
209  }
210  xCoefsAlongY.push_back(yCoefs);
211  }
212  xCoefs.push_back(xCoefsAlongY);
213  }
214  coefficients.push_back(xCoefs);
215  }
216 
217  }
218  else assert(0);
219 
220  RooArgSet leafset;
221  getLeafDependents(theXVar, leafset);
222  getLeafDependents(theYVar, leafset);
223  getLeafDependents(theZVar, leafset);
224  addLeafDependents(leafset);
225 
226  emptyFcnList();
227 }

◆ MELANCSpline_3D_fast() [4/4]

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

Definition at line 229 of file MELANCSpline_3D_fast.cc.

232  :
233  MELANCSplineCore(other, name),
234  rangeYmin(other.rangeYmin), rangeYmax(other.rangeYmax),
235  rangeZmin(other.rangeZmin), rangeZmax(other.rangeZmax),
236  bcBeginX(other.bcBeginX), bcEndX(other.bcEndX),
237  bcBeginY(other.bcBeginY), bcEndY(other.bcEndY),
238  bcBeginZ(other.bcBeginZ), bcEndZ(other.bcEndZ),
239  theYVar("theYVar", this, other.theYVar),
240  theZVar("theZVar", this, other.theZVar),
241  YList(other.YList),
242  ZList(other.ZList),
243  FcnList(other.FcnList),
244  kappaX(other.kappaX),
245  kappaY(other.kappaY),
246  kappaZ(other.kappaZ),
248 {}

◆ ~MELANCSpline_3D_fast()

virtual MELANCSpline_3D_fast::~MELANCSpline_3D_fast ( )
inlinevirtual

Definition at line 69 of file MELANCSpline_3D_fast.h.

69 {}

Member Function Documentation

◆ analyticalIntegral()

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

Implements MELANCSplineCore.

Definition at line 491 of file MELANCSpline_3D_fast.cc.

491  {
492  Double_t value = interpolateFcn(code, rangeName);
493  if (useFloor && value<floorInt){
494  if (verbosity>=MELANCSplineCore::kError) coutE(Integration) << "MELANCSpline_3D_fast ERROR::MELANCSpline_3D_fast(" << GetName() << ") integration returned " << value << " for code = " << code << endl;
495  value = floorInt;
496  }
497  if (verbosity==MELANCSplineCore::kVerbose){ cout << "MELANCSpline_3D_fast(" << GetName() << ")::analyticalIntegral = " << value << " for code = " << code << endl; }
498  return value;
499 }

◆ clone()

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

Implements MELANCSplineCore.

Definition at line 68 of file MELANCSpline_3D_fast.h.

68 { return new MELANCSpline_3D_fast(*this, newname); }

◆ cropValueForRange()

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

Implements MELANCSplineCore.

Definition at line 534 of file MELANCSpline_3D_fast.cc.

534  {
535  if (testRangeValidity(val, whichDirection)) return;
536  const T* range[2];
537  if (whichDirection==0){
538  range[0] = &rangeXmin;
539  range[1] = &rangeXmax;
540  }
541  else if (whichDirection==1){
542  range[0] = &rangeYmin;
543  range[1] = &rangeYmax;
544  }
545  else{
546  range[0] = &rangeZmin;
547  range[1] = &rangeZmax;
548  }
549  if (val<*(range[0])) val = *(range[0]);
550  if (val>*(range[1])) val = *(range[1]);
551 }

◆ emptyFcnList()

virtual void MELANCSpline_3D_fast::emptyFcnList ( )
inlineprotectedvirtual

Implements MELANCSplineCore.

Definition at line 77 of file MELANCSpline_3D_fast.h.

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

◆ evaluate()

Double_t MELANCSpline_3D_fast::evaluate ( ) const
protectedvirtual

Implements MELANCSplineCore.

Definition at line 438 of file MELANCSpline_3D_fast.cc.

438  {
439  Double_t value = interpolateFcn(0);
440  if (useFloor && value<floorEval){
441  if (verbosity>=MELANCSplineCore::kError) coutE(Eval) << "MELANCSpline_3D_fast ERROR::MELANCSpline_3D_fast(" << GetName() << ") evaluation returned " << value << " at (x, y, z) = (" << theXVar << ", " << theYVar << ", " << theZVar << ")" << endl;
442  value = floorEval;
443  }
444  if (verbosity==MELANCSplineCore::kVerbose){ cout << "MELANCSpline_3D_fast(" << GetName() << ")::evaluate = " << value << " at (x, y, z) = (" << theXVar << ", " << theYVar << ", " << theZVar << ")" << endl; }
445  return value;
446 }

◆ getAnalyticalIntegral()

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

Implements MELANCSplineCore.

Definition at line 447 of file MELANCSpline_3D_fast.cc.

447  {
448  if (_forceNumInt) return 0;
449 
450  Int_t code=1;
451 
452  RooArgSet Xdeps, Ydeps, Zdeps;
453  RooRealVar* rrv_x = dynamic_cast<RooRealVar*>(theXVar.absArg());
454  RooRealVar* rrv_y = dynamic_cast<RooRealVar*>(theYVar.absArg());
455  RooRealVar* rrv_z = dynamic_cast<RooRealVar*>(theZVar.absArg());
456  if (rrv_x==0) theXVar.absArg()->leafNodeServerList(&Xdeps, 0, true);
457  if (rrv_y==0) theYVar.absArg()->leafNodeServerList(&Ydeps, 0, true);
458  if (rrv_z==0) theZVar.absArg()->leafNodeServerList(&Zdeps, 0, true);
459 
460  if (rrv_x!=0){
461  if (
462  (Ydeps.find(*rrv_x)==0 || rrv_y!=0)
463  &&
464  (Zdeps.find(*rrv_x)==0 || rrv_z!=0)
465  ){
466  if (matchArgs(allVars, analVars, theXVar)) code*=2;
467  }
468  }
469  if (rrv_y!=0){
470  if (
471  (Xdeps.find(*rrv_y)==0 || rrv_x!=0)
472  &&
473  (Zdeps.find(*rrv_y)==0 || rrv_z!=0)
474  ){
475  if (matchArgs(allVars, analVars, theYVar)) code*=3;
476  }
477  }
478  if (rrv_z!=0){
479  if (
480  (Xdeps.find(*rrv_z)==0 || rrv_x!=0)
481  &&
482  (Ydeps.find(*rrv_z)==0 || rrv_y!=0)
483  ){
484  if (matchArgs(allVars, analVars, theZVar)) code*=5;
485  }
486  }
487 
488  if (code==1) code=0;
489  return code;
490 }

◆ getCoefficientsPerYPerZ()

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

Definition at line 426 of file MELANCSpline_3D_fast.cc.

431  {
432  vector<MELANCSplineCore::T> fcnList;
433  for (unsigned int bin=0; bin<npointsX(); bin++){ fcnList.push_back(FcnList.at(zbin).at(ybin).at(bin)); }
434  vector<vector<MELANCSplineCore::T>> coefs = getCoefficientsAlongDirection(kappaX, xAinv, fcnList, bcBegin, bcEnd, xbin);
435  return coefs;
436 }

◆ getKappas()

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

Implements MELANCSplineCore.

Definition at line 356 of file MELANCSpline_3D_fast.cc.

356  {
357  kappas.clear();
359 
360  Int_t npoints;
361  vector<MELANCSplineCore::T> const* coord;
362  if (whichDirection==0){
363  npoints=npointsX();
364  coord=&XList;
365  }
366  else if (whichDirection==1){
367  npoints=npointsY();
368  coord=&YList;
369  }
370  else{
371  npoints=npointsZ();
372  coord=&ZList;
373  }
374 
375  for (Int_t j=0; j<npoints-1; j++){
376  MELANCSplineCore::T val_j = coord->at(j);
377  MELANCSplineCore::T val_jpo = coord->at(j+1);
378  MELANCSplineCore::T val_diff = (val_jpo-val_j);
379  if (fabs(val_diff)>MELANCSplineCore::T(0)) kappa = 1./val_diff;
380  else kappa = 0;
381  kappas.push_back(kappa);
382  }
383  kappas.push_back(kappa); // Push the same kappa_(N-1)=kappa_(N-2) at the end point
384 }

◆ getTVar()

MELANCSplineCore::T MELANCSpline_3D_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 417 of file MELANCSpline_3D_fast.cc.

417  {
418  const MELANCSplineCore::T& K=kappas.at(bin);
419  vector<MELANCSplineCore::T> const* coord;
420  if (whichDirection==0) coord=&XList;
421  else if (whichDirection==1) coord=&YList;
422  else coord=&ZList;
423  return (val-coord->at(bin))*K;
424 }

◆ getWhichBin()

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

Implements MELANCSplineCore.

Definition at line 385 of file MELANCSpline_3D_fast.cc.

385  {
386  Int_t bin=-1;
387  MELANCSplineCore::T valj, valjpo;
388  Int_t npoints;
389  vector<MELANCSplineCore::T> const* coord;
390  if (whichDirection==0){
391  npoints=npointsX();
392  coord=&XList;
393  }
394  else if (whichDirection==1){
395  npoints=npointsY();
396  coord=&YList;
397  }
398  else{
399  npoints=npointsZ();
400  coord=&ZList;
401  }
402 
403  if (npoints<=1) bin=0;
404  else{
405  valjpo = coord->at(0);
406  for (Int_t j=0; j<npoints-1; j++){
407  valj = coord->at(j);
408  valjpo = coord->at(j+1);
409  if (val<valjpo && val>=valj){ bin=j; break; }
410  }
411  if (bin==-1 && val>=valjpo) bin=npoints-2;
412  else if (bin==-1) bin=0;
413  }
414 
415  return bin;
416 }

◆ interpolateFcn()

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

Implements MELANCSplineCore.

Definition at line 250 of file MELANCSpline_3D_fast.cc.

250  {
252 
253  if (verbosity==MELANCSplineCore::kVerbose) cout << "MELANCSpline_3D_fast(" << GetName() << ")::interpolateFcn begin with code: " << code << endl;
254 
255  // Get bins
256  vector<Int_t> varprime; varprime.push_back(2); varprime.push_back(3); varprime.push_back(5);
257  const Int_t ndims=(const Int_t)varprime.size();
258  vector<Int_t> varbin;
259  vector<Int_t> varbinmin;
260  vector<Int_t> varbinmax;
261  vector<MELANCSplineCore::T> tvar;
262  vector<MELANCSplineCore::T> tvarmin;
263  vector<MELANCSplineCore::T> tvarmax;
264  vector<const vector<MELANCSplineCore::T>*> varkappa; varkappa.push_back(&kappaX); varkappa.push_back(&kappaY); varkappa.push_back(&kappaZ);
265  vector<const RooRealProxy*> varcoord; varcoord.push_back(&theXVar); varcoord.push_back(&theYVar); varcoord.push_back(&theZVar);
266  for (Int_t idim=0; idim<ndims; idim++){
267  Int_t binval=-1;
268  Int_t binmin=-1;
269  Int_t binmax=-1;
270  MELANCSplineCore::T tval=0;
271  MELANCSplineCore::T tmin=0;
272  MELANCSplineCore::T tmax=0;
273  if (code==0 || code%varprime.at(idim)!=0){
274  if (!testRangeValidity(*(varcoord.at(idim)), idim)) return 0;
275  binval = getWhichBin(*(varcoord.at(idim)), idim);
276  tval = getTVar(*(varkappa.at(idim)), *(varcoord.at(idim)), binval, idim);
277  }
278  else{
279  MELANCSplineCore::T coordmin = varcoord.at(idim)->min(rangeName); cropValueForRange(coordmin, idim);
280  MELANCSplineCore::T coordmax = varcoord.at(idim)->max(rangeName); cropValueForRange(coordmax, idim);
281  const std::vector<MELANCSplineCore::T>& kappadim = *(varkappa.at(idim));
282  binmin = getWhichBin(coordmin, idim);
283  tmin = getTVar(kappadim, coordmin, binmin, idim);
284  binmax = getWhichBin(coordmax, idim);
285  tmax = getTVar(kappadim, coordmax, binmax, idim);
286  }
287  varbin.push_back(binval);
288  varbinmin.push_back(binmin);
289  varbinmax.push_back(binmax);
290  tvar.push_back(tval);
291  tvarmin.push_back(tmin);
292  tvarmax.push_back(tmax);
293  }
294 
295  for (int ix=0; ix<(int)coefficients.size(); ix++){
296  if (
297  (varbin[0]>=0 && ix!=varbin[0])
298  ||
299  (varbinmin[0]>=0 && varbinmax[0]>=varbinmin[0] && !(varbinmin[0]<=ix && ix<=varbinmax[0]))
300  ) continue;
301  MELANCSplineCore::T txlow=0, txhigh=1;
302  if (code>0 && code%varprime[0]==0){
303  if (ix==varbinmin[0]) txlow=tvarmin[0];
304  if (ix==varbinmax[0]) txhigh=tvarmax[0];
305  }
306  else txhigh=tvar[0];
307  // Get the x coefficients interpolated across y
308  vector<MELANCSplineCore::T> xCoefs;
309  for (int icx=0; icx<(int)coefficients.at(ix).size(); icx++){
311  for (int iy=0; iy<(int)coefficients.at(ix).at(icx).size(); iy++){
312  if (
313  (varbin[1]>=0 && iy!=varbin[1])
314  ||
315  (varbinmin[1]>=0 && varbinmax[1]>=varbinmin[1] && !(varbinmin[1]<=iy && iy<=varbinmax[1]))
316  ) continue;
317  MELANCSplineCore::T tylow=0, tyhigh=1;
318  if (code>0 && code%varprime[1]==0){
319  if (iy==varbinmin[1]) tylow=tvarmin[1];
320  if (iy==varbinmax[1]) tyhigh=tvarmax[1];
321  }
322  else tyhigh=tvar[1];
323  // Get the y coefficients interpolated across z
324  vector<MELANCSplineCore::T> yCoefs;
325  for (int icy=0; icy<(int)coefficients.at(ix).at(icx).at(iy).size(); icy++){
327  for (int iz=0; iz<(int)coefficients.at(ix).at(icx).at(iy).at(icy).size(); iz++){
328  if (
329  (varbin[2]>=0 && iz!=varbin[2])
330  ||
331  (varbinmin[2]>=0 && varbinmax[2]>=varbinmin[2] && !(varbinmin[2]<=iz && iz<=varbinmax[2]))
332  ) continue;
333 
334  MELANCSplineCore::T tzlow=0, tzhigh=1;
335  if (code>0 && code%varprime[2]==0){
336  if (iz==varbinmin[2]) tzlow=tvarmin[2];
337  if (iz==varbinmax[2]) tzhigh=tvarmax[2];
338  }
339  else tzhigh=tvar[2];
340 
341  theYCoef += evalSplineSegment(coefficients.at(ix).at(icx).at(iy).at(icy).at(iz), varkappa[2]->at(iz), tzhigh, tzlow, (code>0 && code%varprime[2]==0));
342  }
343  yCoefs.push_back(theYCoef);
344  }
345  theXCoef += evalSplineSegment(yCoefs, varkappa[1]->at(iy), tyhigh, tylow, (code>0 && code%varprime[1]==0));
346  }
347  xCoefs.push_back(theXCoef);
348  }
349  // Evaluate value of spline at x with coefficients evaluated at y
350  res += evalSplineSegment(xCoefs, varkappa[0]->at(ix), txhigh, txlow, (code>0 && code%varprime[0]==0));
351  }
352 
353  return res;
354 }

◆ npointsY()

unsigned int MELANCSpline_3D_fast::npointsY ( ) const
inlineprotected

Definition at line 79 of file MELANCSpline_3D_fast.h.

79 { return YList.size(); }

◆ npointsZ()

unsigned int MELANCSpline_3D_fast::npointsZ ( ) const
inlineprotected

Definition at line 80 of file MELANCSpline_3D_fast.h.

80 { return ZList.size(); }

◆ setRangeValidity()

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

Implements MELANCSplineCore.

Definition at line 517 of file MELANCSpline_3D_fast.cc.

517  {
518  T* range[2];
519  if (whichDirection==0){
520  range[0] = &rangeXmin;
521  range[1] = &rangeXmax;
522  }
523  else if (whichDirection==1){
524  range[0] = &rangeYmin;
525  range[1] = &rangeYmax;
526  }
527  else{
528  range[0] = &rangeZmin;
529  range[1] = &rangeZmax;
530  }
531  *(range[0])=valmin;
532  *(range[1])=valmax;
533 }

◆ testRangeValidity()

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

Implements MELANCSplineCore.

Definition at line 501 of file MELANCSpline_3D_fast.cc.

501  {
502  const T* range[2];
503  if (whichDirection==0){
504  range[0] = &rangeXmin;
505  range[1] = &rangeXmax;
506  }
507  else if (whichDirection==1){
508  range[0] = &rangeYmin;
509  range[1] = &rangeYmax;
510  }
511  else{
512  range[0] = &rangeZmin;
513  range[1] = &rangeZmax;
514  }
515  return (*(range[0])>*(range[1]) || (val>=*(range[0]) && val<=*(range[1])));
516 }

Member Data Documentation

◆ bcBeginX

const BoundaryCondition MELANCSpline_3D_fast::bcBeginX
protected

Definition at line 18 of file MELANCSpline_3D_fast.h.

◆ bcBeginY

const BoundaryCondition MELANCSpline_3D_fast::bcBeginY
protected

Definition at line 20 of file MELANCSpline_3D_fast.h.

◆ bcBeginZ

const BoundaryCondition MELANCSpline_3D_fast::bcBeginZ
protected

Definition at line 22 of file MELANCSpline_3D_fast.h.

◆ bcEndX

const BoundaryCondition MELANCSpline_3D_fast::bcEndX
protected

Definition at line 19 of file MELANCSpline_3D_fast.h.

◆ bcEndY

const BoundaryCondition MELANCSpline_3D_fast::bcEndY
protected

Definition at line 21 of file MELANCSpline_3D_fast.h.

◆ bcEndZ

const BoundaryCondition MELANCSpline_3D_fast::bcEndZ
protected

Definition at line 23 of file MELANCSpline_3D_fast.h.

◆ coefficients

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

Definition at line 39 of file MELANCSpline_3D_fast.h.

◆ FcnList

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

Definition at line 30 of file MELANCSpline_3D_fast.h.

◆ kappaX

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

Definition at line 32 of file MELANCSpline_3D_fast.h.

◆ kappaY

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

Definition at line 33 of file MELANCSpline_3D_fast.h.

◆ kappaZ

std::vector<T> MELANCSpline_3D_fast::kappaZ
protected

Definition at line 34 of file MELANCSpline_3D_fast.h.

◆ rangeYmax

T MELANCSpline_3D_fast::rangeYmax
protected

Definition at line 14 of file MELANCSpline_3D_fast.h.

◆ rangeYmin

T MELANCSpline_3D_fast::rangeYmin
protected

Definition at line 13 of file MELANCSpline_3D_fast.h.

◆ rangeZmax

T MELANCSpline_3D_fast::rangeZmax
protected

Definition at line 16 of file MELANCSpline_3D_fast.h.

◆ rangeZmin

T MELANCSpline_3D_fast::rangeZmin
protected

Definition at line 15 of file MELANCSpline_3D_fast.h.

◆ theYVar

RooRealProxy MELANCSpline_3D_fast::theYVar
protected

Definition at line 25 of file MELANCSpline_3D_fast.h.

◆ theZVar

RooRealProxy MELANCSpline_3D_fast::theZVar
protected

Definition at line 26 of file MELANCSpline_3D_fast.h.

◆ YList

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

Definition at line 27 of file MELANCSpline_3D_fast.h.

◆ ZList

std::vector<T> MELANCSpline_3D_fast::ZList
protected

Definition at line 28 of file MELANCSpline_3D_fast.h.


The documentation for this class was generated from the following files:
MELANCSplineCore::kVerbose
@ kVerbose
Definition: MELANCSplineCore.h:25
MELANCSpline_3D_fast::bcEndY
const BoundaryCondition bcEndY
Definition: MELANCSpline_3D_fast.h:21
MELANCSpline_3D_fast::bcBeginZ
const BoundaryCondition bcBeginZ
Definition: MELANCSpline_3D_fast.h:22
MELANCSpline_3D_fast::interpolateFcn
virtual T interpolateFcn(Int_t code, const char *rangeName=0) const
Definition: MELANCSpline_3D_fast.cc:250
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_3D_fast::kappaY
std::vector< T > kappaY
Definition: MELANCSpline_3D_fast.h:33
MELANCSpline_3D_fast::theZVar
RooRealProxy theZVar
Definition: MELANCSpline_3D_fast.h:26
MELANCSpline_3D_fast::FcnList
std::vector< std::vector< std::vector< T > > > FcnList
Definition: MELANCSpline_3D_fast.h:30
MELANCSpline_3D_fast::getKappas
virtual void getKappas(std::vector< T > &kappas, const Int_t whichDirection)
Definition: MELANCSpline_3D_fast.cc:356
MELANCSpline_3D_fast::bcBeginY
const BoundaryCondition bcBeginY
Definition: MELANCSpline_3D_fast.h:20
TNumericUtil::KahanAccumulator
Definition: MELAAccumulators.h:27
MELANCSpline_3D_fast::cropValueForRange
void cropValueForRange(T &val, const Int_t whichDirection) const
Definition: MELANCSpline_3D_fast.cc:534
MELANCSpline_3D_fast::testRangeValidity
Bool_t testRangeValidity(const T &val, const Int_t whichDirection) const
Definition: MELANCSpline_3D_fast.cc:501
MELANCSplineCore::rangeXmax
T rangeXmax
Definition: MELANCSplineCore.h:73
MELANCSpline_3D_fast::theYVar
RooRealProxy theYVar
Definition: MELANCSpline_3D_fast.h:25
MELANCSplineCore::getLeafDependents
void getLeafDependents(RooRealProxy &proxy, RooArgSet &set)
Definition: MELANCSplineCore.cc:275
MELANCSpline_3D_fast::bcBeginX
const BoundaryCondition bcBeginX
Definition: MELANCSpline_3D_fast.h:18
MELANCSplineCore::kError
@ kError
Definition: MELANCSplineCore.h:24
MELANCSpline_3D_fast::getCoefficientsPerYPerZ
virtual std::vector< std::vector< T > > getCoefficientsPerYPerZ(const std::vector< T > &kappaX, const TMatrix_t &xAinv, const Int_t &ybin, const Int_t &zbin, MELANCSplineCore::BoundaryCondition const &bcBegin, MELANCSplineCore::BoundaryCondition const &bcEnd, const Int_t xbin) const
Definition: MELANCSpline_3D_fast.cc:426
MELANCSpline_3D_fast::ZList
std::vector< T > ZList
Definition: MELANCSpline_3D_fast.h:28
MELANCSpline_3D_fast::emptyFcnList
virtual void emptyFcnList()
Definition: MELANCSpline_3D_fast.h:77
MELANCSpline_3D_fast::npointsZ
unsigned int npointsZ() const
Definition: MELANCSpline_3D_fast.h:80
MELANCSpline_3D_fast::MELANCSpline_3D_fast
MELANCSpline_3D_fast()
MELANCSpline_3D_fast::rangeYmin
T rangeYmin
Definition: MELANCSpline_3D_fast.h:13
MELANCSplineCore::rangeXmin
T rangeXmin
Definition: MELANCSplineCore.h:72
MELANCSpline_3D_fast::kappaX
std::vector< T > kappaX
Definition: MELANCSpline_3D_fast.h:32
MELANCSpline_3D_fast::bcEndZ
const BoundaryCondition bcEndZ
Definition: MELANCSpline_3D_fast.h:23
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
testME_all.int
int
Definition: testME_all.py:13
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_3D_fast::rangeZmin
T rangeZmin
Definition: MELANCSpline_3D_fast.h:15
MELANCSplineCore::MELANCSplineCore
MELANCSplineCore()
MELANCSplineCore::npointsX
unsigned int npointsX() const
Definition: MELANCSplineCore.h:84
MELANCSpline_3D_fast::getTVar
virtual T getTVar(const std::vector< T > &kappas, const T &val, const Int_t &bin, const Int_t whichDirection) const
Definition: MELANCSpline_3D_fast.cc:417
MELANCSpline_3D_fast::rangeZmax
T rangeZmax
Definition: MELANCSpline_3D_fast.h:16
dd_global::cout
integer cout
Definition: DD_global.F90:21
ndims
@ ndims
Definition: TMCFM.hh:22
MELANCSplineCore::useFloor
Bool_t useFloor
Definition: MELANCSplineCore.h:68
MELANCSpline_3D_fast::coefficients
std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< T > > > > > > coefficients
Definition: MELANCSpline_3D_fast.h:39
MELANCSpline_3D_fast::npointsY
unsigned int npointsY() const
Definition: MELANCSpline_3D_fast.h:79
MELANCSplineCore::bcNaturalSpline
@ bcNaturalSpline
Definition: MELANCSplineCore.h:32
MELANCSplineCore::floorInt
T floorInt
Definition: MELANCSplineCore.h:70
MELANCSpline_3D_fast::kappaZ
std::vector< T > kappaZ
Definition: MELANCSpline_3D_fast.h:34
MELANCSplineCore::XList
std::vector< T > XList
Definition: MELANCSplineCore.h:77
MELANCSplineCore::floorEval
T floorEval
Definition: MELANCSplineCore.h:69
MELANCSpline_3D_fast::rangeYmax
T rangeYmax
Definition: MELANCSpline_3D_fast.h:14
MELANCSpline_3D_fast::bcEndX
const BoundaryCondition bcEndX
Definition: MELANCSpline_3D_fast.h:19
MELANCSplineCore::addLeafDependents
void addLeafDependents(RooArgSet &set)
Definition: MELANCSplineCore.cc:280
MELANCSpline_3D_fast::getWhichBin
virtual Int_t getWhichBin(const T &val, const Int_t whichDirection) const
Definition: MELANCSpline_3D_fast.cc:385
MELANCSplineCore::verbosity
VerbosityLevel verbosity
Definition: MELANCSplineCore.h:67
MELANCSplineCore::theXVar
RooRealProxy theXVar
Definition: MELANCSplineCore.h:75
MELANCSpline_3D_fast::YList
std::vector< T > YList
Definition: MELANCSpline_3D_fast.h:27
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