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

#include <MELADifermionResolutionModel.h>

Collaboration diagram for MELADifermionResolutionModel:

Public Member Functions

 MELADifermionResolutionModel (TVar::Production prod_, float sqrts_, TString filename, TString extraName)
 
 ~MELADifermionResolutionModel ()
 
float getVal (float val)
 
bool isProduction (TVar::Production prod_)
 
bool isValid ()
 

Protected Attributes

bool valid
 
TVar::Production prod
 
float sqrts
 
RooRealVar * varReco
 
MELANCSplineFactory_1DsplineFactory
 
MELAFuncPdfrecoBW
 

Detailed Description

Definition at line 8 of file MELADifermionResolutionModel.h.

Constructor & Destructor Documentation

◆ MELADifermionResolutionModel()

MELADifermionResolutionModel::MELADifermionResolutionModel ( TVar::Production  prod_,
float  sqrts_,
TString  filename,
TString  extraName 
)

Definition at line 9 of file MELADifermionResolutionModel.cc.

9  :
10 valid(false),
11 prod(prod_),
12 sqrts(sqrts_),
13 varReco(0),
14 splineFactory(0),
15 recoBW(0)
16 {
17  TFile* finput = TFile::Open(filename, "read");
18  if (finput!=0){
19  TGraph* tg = (TGraph*)finput->Get("tg_recoBW");
20  if (tg!=0){
21  TString appendName = extraName + TVar::ProductionName(prod);
22  varReco = new RooRealVar(Form("%s_varReco", appendName.Data()), "", 0.1, sqrts*1000.);
23  splineFactory = new MELANCSplineFactory_1D(*varReco, appendName);
26  valid=(recoBW!=0);
27  }
28  finput->Close();
29  }
30 }

◆ ~MELADifermionResolutionModel()

MELADifermionResolutionModel::~MELADifermionResolutionModel ( )

Definition at line 31 of file MELADifermionResolutionModel.cc.

31  {
32  // Do not delete recoBW. It is owned by splineFactory.
33  delete splineFactory;
34  delete varReco;
35 }

Member Function Documentation

◆ getVal()

float MELADifermionResolutionModel::getVal ( float  val)

Definition at line 36 of file MELADifermionResolutionModel.cc.

36  {
37  float result = -1;
38  if (recoBW!=0){
39  if (val<0.){ /*result=-1;*/ }
40  else if (val<0.1) result=0;
41  else{
42  varReco->setVal(val);
43  result = recoBW->getVal();
44  result /= 2.*val; // Result is BW(s), not BW(m)
45  }
46  }
47  return result;
48 }

◆ isProduction()

bool MELADifermionResolutionModel::isProduction ( TVar::Production  prod_)

Definition at line 49 of file MELADifermionResolutionModel.cc.

49 { return (prod_==prod); }

◆ isValid()

bool MELADifermionResolutionModel::isValid ( )
inline

Definition at line 22 of file MELADifermionResolutionModel.h.

22 { return valid; }

Member Data Documentation

◆ prod

TVar::Production MELADifermionResolutionModel::prod
protected

Definition at line 11 of file MELADifermionResolutionModel.h.

◆ recoBW

MELAFuncPdf* MELADifermionResolutionModel::recoBW
protected

Definition at line 15 of file MELADifermionResolutionModel.h.

◆ splineFactory

MELANCSplineFactory_1D* MELADifermionResolutionModel::splineFactory
protected

Definition at line 14 of file MELADifermionResolutionModel.h.

◆ sqrts

float MELADifermionResolutionModel::sqrts
protected

Definition at line 12 of file MELADifermionResolutionModel.h.

◆ valid

bool MELADifermionResolutionModel::valid
protected

Definition at line 10 of file MELADifermionResolutionModel.h.

◆ varReco

RooRealVar* MELADifermionResolutionModel::varReco
protected

Definition at line 13 of file MELADifermionResolutionModel.h.


The documentation for this class was generated from the following files:
MELADifermionResolutionModel::recoBW
MELAFuncPdf * recoBW
Definition: MELADifermionResolutionModel.h:15
MELADifermionResolutionModel::splineFactory
MELANCSplineFactory_1D * splineFactory
Definition: MELADifermionResolutionModel.h:14
TVar::ProductionName
TString ProductionName(TVar::Production temp)
Definition: TVar.cc:64
MELANCSplineFactory_1D
Definition: MELANCSplineFactory_1D.h:12
MELADifermionResolutionModel::valid
bool valid
Definition: MELADifermionResolutionModel.h:10
MELADifermionResolutionModel::sqrts
float sqrts
Definition: MELADifermionResolutionModel.h:12
MELADifermionResolutionModel::varReco
RooRealVar * varReco
Definition: MELADifermionResolutionModel.h:13
MELANCSplineFactory_1D::setPoints
void setPoints(TTree *tree)
Definition: MELANCSplineFactory_1D.cc:21
MELADifermionResolutionModel::prod
TVar::Production prod
Definition: MELADifermionResolutionModel.h:11
MELANCSplineFactory_1D::getPDF
MELAFuncPdf * getPDF()
Definition: MELANCSplineFactory_1D.h:37