JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
MELANCSplineFactory_1D.cc
Go to the documentation of this file.
2 #include <cassert>
3 
4 using namespace std;
5 
6 
8  RooAbsReal& splineVar_, TString appendName_,
11 ) :
12  appendName(appendName_),
13  bcBeginX(bcBeginX_), bcEndX(bcEndX_),
14  splineVar(&splineVar_),
15  fcn(0),
16  PDF(0)
17 {}
19  destroyPDF();
20 }
22  vector<pair<MELANCSplineCore::T, MELANCSplineCore::T>> pList;
24  tree->SetBranchAddress("X", &x);
25  tree->SetBranchAddress("Fcn", &fcn);
26  int n = tree->GetEntries();
27  for (int ip=0; ip<n; ip++){ tree->GetEntry(ip); pList.push_back(pair<MELANCSplineCore::T, MELANCSplineCore::T>(x, fcn)); }
28  setPoints(pList);
29 }
31  vector<pair<MELANCSplineCore::T, MELANCSplineCore::T>> pList;
32  double* xx = tg->GetX();
33  double* yy = tg->GetY();
34  int n = tg->GetN();
35  for (int ip=0; ip<n; ip++) pList.push_back(pair<MELANCSplineCore::T, MELANCSplineCore::T>(xx[ip], yy[ip]));
36  setPoints(pList);
37 }
38 const std::vector<std::pair<MELANCSplineCore::T, MELANCSplineCore::T>> MELANCSplineFactory_1D::getPoints(const std::vector<MELANCSplineCore::T>& XList, const std::vector<MELANCSplineCore::T>& FcnList){
39  const unsigned int nX = XList.size();
40  const unsigned int n = FcnList.size();
41  if (nX!=n){
42  cerr << "MELANCSplineFactory_1D::getPoints: nX=" << nX << " != nFcn=" << n << endl;
43  assert(0);
44  }
45  std::vector<std::pair<MELANCSplineCore::T, MELANCSplineCore::T>> pList; pList.reserve(n);
46  for (unsigned int ip=0; ip<n; ip++) pList.push_back(pair<MELANCSplineCore::T, MELANCSplineCore::T>(XList.at(ip), FcnList.at(ip)));
47  return pList;
48 }
49 
50 void MELANCSplineFactory_1D::destroyPDF(){ delete PDF; PDF=0; delete fcn; fcn=0; }
51 void MELANCSplineFactory_1D::initPDF(const std::vector<std::pair<MELANCSplineCore::T, MELANCSplineCore::T>>& pList){
52  destroyPDF();
53 
54  const unsigned int n = pList.size();
55  std::vector<MELANCSplineCore::T> XList;
56  std::vector<MELANCSplineCore::T> FcnList;
57  for (unsigned int ip=0; ip<n; ip++){
58  XList.push_back(pList.at(ip).first);
59  FcnList.push_back(pList.at(ip).second);
60  }
61 
62  TString name = "Func";
63  if (appendName!="") name = Form("%s_%s", name.Data(), appendName.Data());
64  TString title=name;
66  name.Data(),
67  title.Data(),
68  *splineVar,
69  XList, FcnList,
71  );
72 
73  name.Prepend("PDF_"); title=name;
74  PDF = new MELAFuncPdf(
75  name.Data(),
76  title.Data(),
77  *fcn
78  );
79 }
80 
84  const unsigned int /*direction*/
85 ){
86  bcBeginX=bcBegin;
87  bcEndX=bcEnd;
88 }
MELANCSplineFactory_1D::MELANCSplineFactory_1D
MELANCSplineFactory_1D(RooAbsReal &splineVar_, TString appendName_="", MELANCSplineCore::BoundaryCondition const bcBeginX_=MELANCSplineCore::bcNaturalSpline, MELANCSplineCore::BoundaryCondition const bcEndX_=MELANCSplineCore::bcNaturalSpline)
Definition: MELANCSplineFactory_1D.cc:7
MELANCSplineFactory_1D::setEndConditions
void setEndConditions(MELANCSplineCore::BoundaryCondition const bcBegin, MELANCSplineCore::BoundaryCondition const bcEnd, const unsigned int direction=0)
Definition: MELANCSplineFactory_1D.cc:81
MELANCSplineCore::T
Float_t T
Definition: MELANCSplineCore.h:18
MELANCSplineFactory_1D::initPDF
void initPDF(const std::vector< std::pair< MELANCSplineCore::T, MELANCSplineCore::T >> &pList)
Definition: MELANCSplineFactory_1D.cc:51
MELANCSplineFactory_1D.h
MELANCSplineFactory_1D::splineVar
RooAbsReal * splineVar
Definition: MELANCSplineFactory_1D.h:19
MELANCSplineFactory_1D::bcEndX
MELANCSplineCore::BoundaryCondition bcEndX
Definition: MELANCSplineFactory_1D.h:17
MELANCSplineCore::BoundaryCondition
BoundaryCondition
Definition: MELANCSplineCore.h:28
MELANCSplineFactory_1D::bcBeginX
MELANCSplineCore::BoundaryCondition bcBeginX
Definition: MELANCSplineFactory_1D.h:16
MELANCSplineFactory_1D::fcn
MELANCSpline_1D_fast * fcn
Definition: MELANCSplineFactory_1D.h:20
MELANCSpline_1D_fast
Definition: MELANCSpline_1D_fast.h:12
MELANCSplineFactory_1D::getPoints
const std::vector< std::pair< MELANCSplineCore::T, MELANCSplineCore::T > > getPoints(const std::vector< MELANCSplineCore::T > &XList, const std::vector< MELANCSplineCore::T > &FcnList)
Definition: MELANCSplineFactory_1D.cc:38
MELAFuncPdf
Definition: MELAFuncPdf.h:7
MELANCSplineFactory_1D::PDF
MELAFuncPdf * PDF
Definition: MELANCSplineFactory_1D.h:21
MELANCSplineFactory_1D::setPoints
void setPoints(TTree *tree)
Definition: MELANCSplineFactory_1D.cc:21
MELANCSplineFactory_1D::destroyPDF
void destroyPDF()
Definition: MELANCSplineFactory_1D.cc:50
MELANCSplineFactory_1D::~MELANCSplineFactory_1D
~MELANCSplineFactory_1D()
Definition: MELANCSplineFactory_1D.cc:18
MELANCSplineFactory_1D::appendName
TString appendName
Definition: MELANCSplineFactory_1D.h:14