JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
MELANCSplineFactory_3D.h
Go to the documentation of this file.
1 #ifndef MELANCSPLINEFACTORY_3D
2 #define MELANCSPLINEFACTORY_3D
3 
4 #include <vector>
5 #include <utility>
6 #include <algorithm>
7 #include "TTree.h"
8 #include "TNumericUtil.hh"
9 #include "MELANCSpline_3D_fast.h"
10 #include "MELAFuncPdf.h"
11 
13 
15 protected:
16  TString appendName;
17 
24 
25  RooAbsReal* XVar;
26  RooAbsReal* YVar;
27  RooAbsReal* ZVar;
30 
31  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);
32 
33  void destroyPDF();
34  void initPDF(const std::vector<splineQuadruplet_t>& pList);
35 
36  void addUnique(std::vector<MELANCSplineCore::T>& list, MELANCSplineCore::T val);
37 
38 public:
40  RooAbsReal& XVar_, RooAbsReal& YVar_, RooAbsReal& ZVar_, TString appendName_="",
47  );
49 
51  MELAFuncPdf* getPDF(){ return PDF; }
52 
53  void setEndConditions(
56  const unsigned int direction
57  );
58 
59  void setPoints(TTree* tree);
60  void setPoints(const std::vector<splineQuadruplet_t>& pList){ initPDF(pList); }
61  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){
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  }
73 
74 };
75 
76 template void MELANCSplineFactory_3D::setPoints<Float_t>(const std::vector<Float_t>& XList, const std::vector<Float_t>& YList, const std::vector<Float_t>& ZList, const std::vector<Float_t>& FcnList);
77 template void MELANCSplineFactory_3D::setPoints<Double_t>(const std::vector<Double_t>& XList, const std::vector<Double_t>& YList, const std::vector<Double_t>& ZList, const std::vector<Double_t>& FcnList);
78 
79 
80 #endif
81 
82 
83 
MELANCSplineFactory_3D::getPDF
MELAFuncPdf * getPDF()
Definition: MELANCSplineFactory_3D.h:51
MELANCSplineFactory_3D::setPoints
void setPoints(const std::vector< splineQuadruplet_t > &pList)
Definition: MELANCSplineFactory_3D.h:60
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
TNumericUtil.hh
MELANCSplineFactory_3D::addUnique
void addUnique(std::vector< MELANCSplineCore::T > &list, MELANCSplineCore::T val)
Definition: MELANCSplineFactory_3D.cc:27
MELANCSplineCore::BoundaryCondition
BoundaryCondition
Definition: MELANCSplineCore.h:28
MELANCSplineFactory_3D::fcn
MELANCSpline_3D_fast * fcn
Definition: MELANCSplineFactory_3D.h:28
MELANCSplineFactory_3D::~MELANCSplineFactory_3D
~MELANCSplineFactory_3D()
Definition: MELANCSplineFactory_3D.cc:23
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
Definition: MELANCSplineFactory_3D.h:14
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::setEndConditions
void setEndConditions(MELANCSplineCore::BoundaryCondition const bcBegin, MELANCSplineCore::BoundaryCondition const bcEnd, const unsigned int direction)
Definition: MELANCSplineFactory_3D.cc:124
MELANCSplineFactory_3D::bcEndX
MELANCSplineCore::BoundaryCondition bcEndX
Definition: MELANCSplineFactory_3D.h:19
MELANCSplineFactory_3D::bcBeginZ
MELANCSplineCore::BoundaryCondition bcBeginZ
Definition: MELANCSplineFactory_3D.h:22
MELAFuncPdf.h
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
MELANCSplineFactory_3D::setPoints
void setPoints(const std::vector< inType > &XList, const std::vector< inType > &YList, const std::vector< inType > &ZList, const std::vector< inType > &FcnList)
Definition: MELANCSplineFactory_3D.h:61
MELANCSplineCore::bcNaturalSpline
@ bcNaturalSpline
Definition: MELANCSplineCore.h:32
MELANCSplineFactory_3D::getFunc
MELANCSpline_3D_fast * getFunc()
Definition: MELANCSplineFactory_3D.h:50
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
TNumericUtil::quadruplet
Definition: TNumericUtil.hh:27
MELANCSplineFactory_3D::setPoints
void setPoints(TTree *tree)
Definition: MELANCSplineFactory_3D.cc:112
MELANCSpline_3D_fast.h
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: MELANCSplineFactory_3D.cc:6