JHUGen MELA
JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
MELADifermionResolutionModel.cc
Go to the documentation of this file.
1
#include <fstream>
2
#include <string>
3
#include <vector>
4
#include <cassert>
5
#include "TGraph.h"
6
#include "
MELADifermionResolutionModel.h
"
7
8
9
MELADifermionResolutionModel::MELADifermionResolutionModel
(
TVar::Production
prod_,
float
sqrts_, TString filename, TString extraName) :
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);
24
splineFactory
->
setPoints
(tg);
25
recoBW
=
splineFactory
->
getPDF
();
26
valid
=(
recoBW
!=0);
27
}
28
finput->Close();
29
}
30
}
31
MELADifermionResolutionModel::~MELADifermionResolutionModel
(){
32
// Do not delete recoBW. It is owned by splineFactory.
33
delete
splineFactory
;
34
delete
varReco
;
35
}
36
float
MELADifermionResolutionModel::getVal
(
float
val){
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
}
49
bool
MELADifermionResolutionModel::isProduction
(
TVar::Production
prod_){
return
(prod_==
prod
); }
50
MELADifermionResolutionModel.h
MELADifermionResolutionModel::getVal
float getVal(float val)
Definition:
MELADifermionResolutionModel.cc:36
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::isProduction
bool isProduction(TVar::Production prod_)
Definition:
MELADifermionResolutionModel.cc:49
MELADifermionResolutionModel::sqrts
float sqrts
Definition:
MELADifermionResolutionModel.h:12
MELADifermionResolutionModel::MELADifermionResolutionModel
MELADifermionResolutionModel(TVar::Production prod_, float sqrts_, TString filename, TString extraName)
Definition:
MELADifermionResolutionModel.cc:9
MELADifermionResolutionModel::varReco
RooRealVar * varReco
Definition:
MELADifermionResolutionModel.h:13
MELANCSplineFactory_1D::setPoints
void setPoints(TTree *tree)
Definition:
MELANCSplineFactory_1D.cc:21
MELADifermionResolutionModel::~MELADifermionResolutionModel
~MELADifermionResolutionModel()
Definition:
MELADifermionResolutionModel.cc:31
MELADifermionResolutionModel::prod
TVar::Production prod
Definition:
MELADifermionResolutionModel.h:11
sqrts
double sqrts
Definition:
TMCFM.hh:290
TVar::Production
Production
Definition:
TVar.hh:61
MELANCSplineFactory_1D::getPDF
MELAFuncPdf * getPDF()
Definition:
MELANCSplineFactory_1D.h:37
MELA
src
MELADifermionResolutionModel.cc
Generated on Fri Oct 18 2024 15:51:14 for JHUGen MELA by
1.8.17