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

#include <MELANCSplineFactory_1D.h>

Collaboration diagram for MELANCSplineFactory_1D:

Public Member Functions

 MELANCSplineFactory_1D (RooAbsReal &splineVar_, TString appendName_="", MELANCSplineCore::BoundaryCondition const bcBeginX_=MELANCSplineCore::bcNaturalSpline, MELANCSplineCore::BoundaryCondition const bcEndX_=MELANCSplineCore::bcNaturalSpline)
 
 ~MELANCSplineFactory_1D ()
 
MELANCSpline_1D_fastgetFunc ()
 
MELAFuncPdfgetPDF ()
 
void setEndConditions (MELANCSplineCore::BoundaryCondition const bcBegin, MELANCSplineCore::BoundaryCondition const bcEnd, const unsigned int direction=0)
 
void setPoints (TTree *tree)
 
void setPoints (TGraph *tg)
 
void setPoints (const std::vector< std::pair< MELANCSplineCore::T, MELANCSplineCore::T >> &pList)
 
template<typename inType >
void setPoints (const std::vector< inType > &XList, const std::vector< inType > &FcnList)
 

Protected Member Functions

const std::vector< std::pair< MELANCSplineCore::T, MELANCSplineCore::T > > getPoints (const std::vector< MELANCSplineCore::T > &XList, const std::vector< MELANCSplineCore::T > &FcnList)
 
void destroyPDF ()
 
void initPDF (const std::vector< std::pair< MELANCSplineCore::T, MELANCSplineCore::T >> &pList)
 

Protected Attributes

TString appendName
 
MELANCSplineCore::BoundaryCondition bcBeginX
 
MELANCSplineCore::BoundaryCondition bcEndX
 
RooAbsReal * splineVar
 
MELANCSpline_1D_fastfcn
 
MELAFuncPdfPDF
 

Detailed Description

Definition at line 12 of file MELANCSplineFactory_1D.h.

Constructor & Destructor Documentation

◆ MELANCSplineFactory_1D()

MELANCSplineFactory_1D::MELANCSplineFactory_1D ( RooAbsReal &  splineVar_,
TString  appendName_ = "",
MELANCSplineCore::BoundaryCondition const  bcBeginX_ = MELANCSplineCore::bcNaturalSpline,
MELANCSplineCore::BoundaryCondition const  bcEndX_ = MELANCSplineCore::bcNaturalSpline 
)

Definition at line 7 of file MELANCSplineFactory_1D.cc.

11  :
12  appendName(appendName_),
13  bcBeginX(bcBeginX_), bcEndX(bcEndX_),
14  splineVar(&splineVar_),
15  fcn(0),
16  PDF(0)
17 {}

◆ ~MELANCSplineFactory_1D()

MELANCSplineFactory_1D::~MELANCSplineFactory_1D ( )

Definition at line 18 of file MELANCSplineFactory_1D.cc.

18  {
19  destroyPDF();
20 }

Member Function Documentation

◆ destroyPDF()

void MELANCSplineFactory_1D::destroyPDF ( )
protected

Definition at line 50 of file MELANCSplineFactory_1D.cc.

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

◆ getFunc()

MELANCSpline_1D_fast* MELANCSplineFactory_1D::getFunc ( )
inline

Definition at line 36 of file MELANCSplineFactory_1D.h.

36 { return fcn; }

◆ getPDF()

MELAFuncPdf* MELANCSplineFactory_1D::getPDF ( )
inline

Definition at line 37 of file MELANCSplineFactory_1D.h.

37 { return PDF; }

◆ getPoints()

const std::vector< std::pair< MELANCSplineCore::T, MELANCSplineCore::T > > MELANCSplineFactory_1D::getPoints ( const std::vector< MELANCSplineCore::T > &  XList,
const std::vector< MELANCSplineCore::T > &  FcnList 
)
protected

Definition at line 38 of file MELANCSplineFactory_1D.cc.

38  {
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 }

◆ initPDF()

void MELANCSplineFactory_1D::initPDF ( const std::vector< std::pair< MELANCSplineCore::T, MELANCSplineCore::T >> &  pList)
protected

Definition at line 51 of file MELANCSplineFactory_1D.cc.

51  {
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 }

◆ setEndConditions()

void MELANCSplineFactory_1D::setEndConditions ( MELANCSplineCore::BoundaryCondition const  bcBegin,
MELANCSplineCore::BoundaryCondition const  bcEnd,
const unsigned int  direction = 0 
)

Definition at line 81 of file MELANCSplineFactory_1D.cc.

85  {
86  bcBeginX=bcBegin;
87  bcEndX=bcEnd;
88 }

◆ setPoints() [1/4]

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

Definition at line 48 of file MELANCSplineFactory_1D.h.

48  {
49  std::vector<MELANCSplineCore::T> transXList;
50  std::vector<MELANCSplineCore::T> transFcnList;
51  for (unsigned int ip=0; ip<XList.size(); ip++) transXList.push_back((MELANCSplineCore::T)XList.at(ip));
52  for (unsigned int ip=0; ip<FcnList.size(); ip++) transFcnList.push_back((MELANCSplineCore::T)FcnList.at(ip));
53  const std::vector<std::pair<MELANCSplineCore::T, MELANCSplineCore::T>> pList = getPoints(transXList, transFcnList);
54  initPDF(pList);
55  }

◆ setPoints() [2/4]

void MELANCSplineFactory_1D::setPoints ( const std::vector< std::pair< MELANCSplineCore::T, MELANCSplineCore::T >> &  pList)
inline

Definition at line 47 of file MELANCSplineFactory_1D.h.

47 { initPDF(pList); }

◆ setPoints() [3/4]

void MELANCSplineFactory_1D::setPoints ( TGraph *  tg)

Definition at line 30 of file MELANCSplineFactory_1D.cc.

30  {
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 }

◆ setPoints() [4/4]

void MELANCSplineFactory_1D::setPoints ( TTree *  tree)

Definition at line 21 of file MELANCSplineFactory_1D.cc.

21  {
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 }

Member Data Documentation

◆ appendName

TString MELANCSplineFactory_1D::appendName
protected

Definition at line 14 of file MELANCSplineFactory_1D.h.

◆ bcBeginX

MELANCSplineCore::BoundaryCondition MELANCSplineFactory_1D::bcBeginX
protected

Definition at line 16 of file MELANCSplineFactory_1D.h.

◆ bcEndX

MELANCSplineCore::BoundaryCondition MELANCSplineFactory_1D::bcEndX
protected

Definition at line 17 of file MELANCSplineFactory_1D.h.

◆ fcn

MELANCSpline_1D_fast* MELANCSplineFactory_1D::fcn
protected

Definition at line 20 of file MELANCSplineFactory_1D.h.

◆ PDF

MELAFuncPdf* MELANCSplineFactory_1D::PDF
protected

Definition at line 21 of file MELANCSplineFactory_1D.h.

◆ splineVar

RooAbsReal* MELANCSplineFactory_1D::splineVar
protected

Definition at line 19 of file MELANCSplineFactory_1D.h.


The documentation for this class was generated from the following files:
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::splineVar
RooAbsReal * splineVar
Definition: MELANCSplineFactory_1D.h:19
MELANCSplineFactory_1D::bcEndX
MELANCSplineCore::BoundaryCondition bcEndX
Definition: MELANCSplineFactory_1D.h:17
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::appendName
TString appendName
Definition: MELANCSplineFactory_1D.h:14