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

#include <MELANCSplineFactory_3D.h>

Collaboration diagram for MELANCSplineFactory_3D:

Public Member Functions

 MELANCSplineFactory_3D (RooAbsReal &XVar_, RooAbsReal &YVar_, RooAbsReal &ZVar_, TString appendName_="", 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)
 
 ~MELANCSplineFactory_3D ()
 
MELANCSpline_3D_fastgetFunc ()
 
MELAFuncPdfgetPDF ()
 
void setEndConditions (MELANCSplineCore::BoundaryCondition const bcBegin, MELANCSplineCore::BoundaryCondition const bcEnd, const unsigned int direction)
 
void setPoints (TTree *tree)
 
void setPoints (const std::vector< splineQuadruplet_t > &pList)
 
template<typename inType >
void setPoints (const std::vector< inType > &XList, const std::vector< inType > &YList, const std::vector< inType > &ZList, const std::vector< inType > &FcnList)
 

Protected Member Functions

const std::vector< splineQuadruplet_tgetPoints (const std::vector< MELANCSplineCore::T > &XList, const std::vector< MELANCSplineCore::T > &YList, const std::vector< MELANCSplineCore::T > &ZList, const std::vector< MELANCSplineCore::T > &FcnList)
 
void destroyPDF ()
 
void initPDF (const std::vector< splineQuadruplet_t > &pList)
 
void addUnique (std::vector< MELANCSplineCore::T > &list, MELANCSplineCore::T val)
 

Protected Attributes

TString appendName
 
MELANCSplineCore::BoundaryCondition bcBeginX
 
MELANCSplineCore::BoundaryCondition bcEndX
 
MELANCSplineCore::BoundaryCondition bcBeginY
 
MELANCSplineCore::BoundaryCondition bcEndY
 
MELANCSplineCore::BoundaryCondition bcBeginZ
 
MELANCSplineCore::BoundaryCondition bcEndZ
 
RooAbsReal * XVar
 
RooAbsReal * YVar
 
RooAbsReal * ZVar
 
MELANCSpline_3D_fastfcn
 
MELAFuncPdfPDF
 

Detailed Description

Definition at line 14 of file MELANCSplineFactory_3D.h.

Constructor & Destructor Documentation

◆ MELANCSplineFactory_3D()

MELANCSplineFactory_3D::MELANCSplineFactory_3D ( RooAbsReal &  XVar_,
RooAbsReal &  YVar_,
RooAbsReal &  ZVar_,
TString  appendName_ = "",
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 
)

Definition at line 6 of file MELANCSplineFactory_3D.cc.

14  :
15  appendName(appendName_),
16  bcBeginX(bcBeginX_), bcEndX(bcEndX_),
17  bcBeginY(bcBeginY_), bcEndY(bcEndY_),
18  bcBeginZ(bcBeginZ_), bcEndZ(bcEndZ_),
19  XVar(&XVar_), YVar(&YVar_), ZVar(&ZVar_),
20  fcn(0),
21  PDF(0)
22 {}

◆ ~MELANCSplineFactory_3D()

MELANCSplineFactory_3D::~MELANCSplineFactory_3D ( )

Definition at line 23 of file MELANCSplineFactory_3D.cc.

23  {
24  destroyPDF();
25 }

Member Function Documentation

◆ addUnique()

void MELANCSplineFactory_3D::addUnique ( std::vector< MELANCSplineCore::T > &  list,
MELANCSplineCore::T  val 
)
protected

Definition at line 27 of file MELANCSplineFactory_3D.cc.

27  {
28  for (unsigned int ip=0; ip<list.size(); ip++){ if (list.at(ip)==val) return; }
29  list.push_back(val);
30 }

◆ destroyPDF()

void MELANCSplineFactory_3D::destroyPDF ( )
protected

Definition at line 61 of file MELANCSplineFactory_3D.cc.

61 { delete PDF; PDF=0; delete fcn; fcn=0; }

◆ getFunc()

MELANCSpline_3D_fast* MELANCSplineFactory_3D::getFunc ( )
inline

Definition at line 50 of file MELANCSplineFactory_3D.h.

50 { return fcn; }

◆ getPDF()

MELAFuncPdf* MELANCSplineFactory_3D::getPDF ( )
inline

Definition at line 51 of file MELANCSplineFactory_3D.h.

51 { return PDF; }

◆ getPoints()

const std::vector< splineQuadruplet_t > MELANCSplineFactory_3D::getPoints ( const std::vector< MELANCSplineCore::T > &  XList,
const std::vector< MELANCSplineCore::T > &  YList,
const std::vector< MELANCSplineCore::T > &  ZList,
const std::vector< MELANCSplineCore::T > &  FcnList 
)
protected

Definition at line 31 of file MELANCSplineFactory_3D.cc.

36  {
37  const unsigned int nX = XList.size();
38  const unsigned int nY = YList.size();
39  const unsigned int nZ = ZList.size();
40  const unsigned int n = FcnList.size();
41  if (nX*nY*nZ!=n){
42  cerr << "MELANCSplineFactory_3D::getPoints: nX=" << nX << " x nY=" << nY << " x nZ=" << nZ << " != nFcn=" << n << endl;
43  assert(0);
44  }
45 
46  std::vector<splineQuadruplet_t> pList; pList.reserve(n);
47  for (unsigned int ix=0; ix<nX; ix++){
48  MELANCSplineCore::T xval = XList.at(ix);
49  for (unsigned int iy=0; iy<nY; iy++){
50  MELANCSplineCore::T yval = YList.at(iy);
51  for (unsigned int iz=0; iz<nZ; iz++){
52  MELANCSplineCore::T zval = ZList.at(iz);
53  unsigned int ip = nZ*(nY*ix + iy) + iz;
54  pList.push_back(splineQuadruplet_t(xval, yval, zval, FcnList.at(ip)));
55  }
56  }
57  }
58  return pList;
59 }

◆ initPDF()

void MELANCSplineFactory_3D::initPDF ( const std::vector< splineQuadruplet_t > &  pList)
protected

Definition at line 62 of file MELANCSplineFactory_3D.cc.

62  {
63  destroyPDF();
64 
65  const unsigned int n = pList.size();
66  vector<MELANCSplineCore::T> XList;
67  vector<MELANCSplineCore::T> YList;
68  vector<MELANCSplineCore::T> ZList;
69  vector<vector<vector<MELANCSplineCore::T>>> FcnList;
70  for (unsigned int ip=0; ip<n; ip++){
71  addUnique(XList, (pList.at(ip))[0]);
72  addUnique(YList, (pList.at(ip))[1]);
73  addUnique(ZList, (pList.at(ip))[2]);
74  }
75  FcnList.reserve(ZList.size());
76  for (unsigned int iz=0; iz<ZList.size(); iz++){
77  vector<vector<MELANCSplineCore::T>> dumz;
78  dumz.reserve(YList.size());
79  for (unsigned int iy=0; iy<YList.size(); iy++){
80  vector<MELANCSplineCore::T> dumy;
81  dumy.reserve(XList.size());
82  for (unsigned int ix=0; ix<XList.size(); ix++){
83  unsigned int ip = ZList.size()*(YList.size()*ix + iy) + iz;
84  dumy.push_back((pList.at(ip))[3]); // Do not use unique here
85  }
86  dumz.push_back(dumy);
87  }
88  FcnList.push_back(dumz);
89  }
90 
91  TString name = "Func";
92  if (appendName!="") name = Form("%s_%s", name.Data(), appendName.Data());
93  TString title=name;
95  name.Data(),
96  title.Data(),
97  *XVar, *YVar, *ZVar,
98  XList, YList, ZList, FcnList,
100  bcBeginY, bcEndY,
102  );
103 
104  name.Prepend("PDF_"); title=name;
105  PDF = new MELAFuncPdf(
106  name.Data(),
107  title.Data(),
108  *fcn
109  );
110 }

◆ setEndConditions()

void MELANCSplineFactory_3D::setEndConditions ( MELANCSplineCore::BoundaryCondition const  bcBegin,
MELANCSplineCore::BoundaryCondition const  bcEnd,
const unsigned int  direction 
)

Definition at line 124 of file MELANCSplineFactory_3D.cc.

128  {
129  switch (direction){
130  case 0:
131  bcBeginX=bcBegin;
132  bcEndX=bcEnd;
133  break;
134  case 1:
135  bcBeginY=bcBegin;
136  bcEndY=bcEnd;
137  break;
138  case 2:
139  bcBeginZ=bcBegin;
140  bcEndZ=bcEnd;
141  break;
142  default:
143  // Do nothing
144  break;
145  }
146 }

◆ setPoints() [1/3]

template<typename inType >
void MELANCSplineFactory_3D::setPoints ( const std::vector< inType > &  XList,
const std::vector< inType > &  YList,
const std::vector< inType > &  ZList,
const std::vector< inType > &  FcnList 
)
inline

Definition at line 61 of file MELANCSplineFactory_3D.h.

61  {
62  std::vector<MELANCSplineCore::T> transXList;
63  std::vector<MELANCSplineCore::T> transYList;
64  std::vector<MELANCSplineCore::T> transZList;
65  std::vector<MELANCSplineCore::T> transFcnList;
66  for (unsigned int ip=0; ip<XList.size(); ip++) transXList.push_back((MELANCSplineCore::T)XList.at(ip));
67  for (unsigned int ip=0; ip<YList.size(); ip++) transYList.push_back((MELANCSplineCore::T)YList.at(ip));
68  for (unsigned int ip=0; ip<ZList.size(); ip++) transZList.push_back((MELANCSplineCore::T)ZList.at(ip));
69  for (unsigned int ip=0; ip<FcnList.size(); ip++) transFcnList.push_back((MELANCSplineCore::T)FcnList.at(ip));
70  const std::vector<splineQuadruplet_t> pList = getPoints(transXList, transYList, transZList, transFcnList);
71  setPoints(pList);
72  }

◆ setPoints() [2/3]

void MELANCSplineFactory_3D::setPoints ( const std::vector< splineQuadruplet_t > &  pList)
inline

Definition at line 60 of file MELANCSplineFactory_3D.h.

60 { initPDF(pList); }

◆ setPoints() [3/3]

void MELANCSplineFactory_3D::setPoints ( TTree *  tree)

Definition at line 112 of file MELANCSplineFactory_3D.cc.

112  {
113  vector<splineQuadruplet_t> pList;
114  MELANCSplineCore::T x, y, z, fcn;
115  tree->SetBranchAddress("X", &x);
116  tree->SetBranchAddress("Y", &y);
117  tree->SetBranchAddress("Z", &z);
118  tree->SetBranchAddress("Fcn", &fcn);
119  int n = tree->GetEntries();
120  for (int ip=0; ip<n; ip++){ tree->GetEntry(ip); pList.push_back(splineQuadruplet_t(x, y, z, fcn)); }
121  setPoints(pList);
122 }

Member Data Documentation

◆ appendName

TString MELANCSplineFactory_3D::appendName
protected

Definition at line 16 of file MELANCSplineFactory_3D.h.

◆ bcBeginX

MELANCSplineCore::BoundaryCondition MELANCSplineFactory_3D::bcBeginX
protected

Definition at line 18 of file MELANCSplineFactory_3D.h.

◆ bcBeginY

MELANCSplineCore::BoundaryCondition MELANCSplineFactory_3D::bcBeginY
protected

Definition at line 20 of file MELANCSplineFactory_3D.h.

◆ bcBeginZ

MELANCSplineCore::BoundaryCondition MELANCSplineFactory_3D::bcBeginZ
protected

Definition at line 22 of file MELANCSplineFactory_3D.h.

◆ bcEndX

MELANCSplineCore::BoundaryCondition MELANCSplineFactory_3D::bcEndX
protected

Definition at line 19 of file MELANCSplineFactory_3D.h.

◆ bcEndY

MELANCSplineCore::BoundaryCondition MELANCSplineFactory_3D::bcEndY
protected

Definition at line 21 of file MELANCSplineFactory_3D.h.

◆ bcEndZ

MELANCSplineCore::BoundaryCondition MELANCSplineFactory_3D::bcEndZ
protected

Definition at line 23 of file MELANCSplineFactory_3D.h.

◆ fcn

MELANCSpline_3D_fast* MELANCSplineFactory_3D::fcn
protected

Definition at line 28 of file MELANCSplineFactory_3D.h.

◆ PDF

MELAFuncPdf* MELANCSplineFactory_3D::PDF
protected

Definition at line 29 of file MELANCSplineFactory_3D.h.

◆ XVar

RooAbsReal* MELANCSplineFactory_3D::XVar
protected

Definition at line 25 of file MELANCSplineFactory_3D.h.

◆ YVar

RooAbsReal* MELANCSplineFactory_3D::YVar
protected

Definition at line 26 of file MELANCSplineFactory_3D.h.

◆ ZVar

RooAbsReal* MELANCSplineFactory_3D::ZVar
protected

Definition at line 27 of file MELANCSplineFactory_3D.h.


The documentation for this class was generated from the following files:
MELANCSplineCore::T
Float_t T
Definition: MELANCSplineCore.h:18
MELANCSplineFactory_3D::appendName
TString appendName
Definition: MELANCSplineFactory_3D.h:16
MELANCSplineFactory_3D::bcEndZ
MELANCSplineCore::BoundaryCondition bcEndZ
Definition: MELANCSplineFactory_3D.h:23
MELANCSplineFactory_3D::bcBeginY
MELANCSplineCore::BoundaryCondition bcBeginY
Definition: MELANCSplineFactory_3D.h:20
MELANCSplineFactory_3D::initPDF
void initPDF(const std::vector< splineQuadruplet_t > &pList)
Definition: MELANCSplineFactory_3D.cc:62
MELANCSplineFactory_3D::addUnique
void addUnique(std::vector< MELANCSplineCore::T > &list, MELANCSplineCore::T val)
Definition: MELANCSplineFactory_3D.cc:27
MELANCSplineFactory_3D::fcn
MELANCSpline_3D_fast * fcn
Definition: MELANCSplineFactory_3D.h:28
MELANCSplineFactory_3D::getPoints
const std::vector< splineQuadruplet_t > getPoints(const std::vector< MELANCSplineCore::T > &XList, const std::vector< MELANCSplineCore::T > &YList, const std::vector< MELANCSplineCore::T > &ZList, const std::vector< MELANCSplineCore::T > &FcnList)
Definition: MELANCSplineFactory_3D.cc:31
MELANCSpline_3D_fast
Definition: MELANCSpline_3D_fast.h:11
MELANCSplineFactory_3D::bcBeginX
MELANCSplineCore::BoundaryCondition bcBeginX
Definition: MELANCSplineFactory_3D.h:18
splineQuadruplet_t
TNumericUtil::quadruplet< MELANCSplineCore::T > splineQuadruplet_t
Definition: MELANCSplineFactory_3D.h:12
MELANCSplineFactory_3D::bcEndX
MELANCSplineCore::BoundaryCondition bcEndX
Definition: MELANCSplineFactory_3D.h:19
MELANCSplineFactory_3D::bcBeginZ
MELANCSplineCore::BoundaryCondition bcBeginZ
Definition: MELANCSplineFactory_3D.h:22
MELANCSplineFactory_3D::YVar
RooAbsReal * YVar
Definition: MELANCSplineFactory_3D.h:26
MELANCSplineFactory_3D::destroyPDF
void destroyPDF()
Definition: MELANCSplineFactory_3D.cc:61
MELANCSplineFactory_3D::ZVar
RooAbsReal * ZVar
Definition: MELANCSplineFactory_3D.h:27
MELANCSplineFactory_3D::bcEndY
MELANCSplineCore::BoundaryCondition bcEndY
Definition: MELANCSplineFactory_3D.h:21
MELAFuncPdf
Definition: MELAFuncPdf.h:7
MELANCSplineFactory_3D::XVar
RooAbsReal * XVar
Definition: MELANCSplineFactory_3D.h:25
MELANCSplineFactory_3D::PDF
MELAFuncPdf * PDF
Definition: MELANCSplineFactory_3D.h:29
globalc::z
double complex, dimension(2, 2) z
Definition: reductionC.F90:50
MELANCSplineFactory_3D::setPoints
void setPoints(TTree *tree)
Definition: MELANCSplineFactory_3D.cc:112