1 |
#include "TROOT.h"
|
2 |
#include "TH1D.h"
|
3 |
#include "RooDataSet.h"
|
4 |
#include "RooRealVar.h"
|
5 |
#include "RooAbsPdf.h"
|
6 |
#include "RooBreitWigner.h"
|
7 |
#include "RooCBShape.h"
|
8 |
#include "RooGaussian.h"
|
9 |
#include "RooFFTConvPdf.h"
|
10 |
#include "RooDataHist.h"
|
11 |
#include "RooHistPdf.h"
|
12 |
#include "RooVoigtianShape.h"
|
13 |
#include "RooKeysPdf.h"
|
14 |
|
15 |
class CSignalModel
|
16 |
{
|
17 |
public:
|
18 |
CSignalModel():model(0){}
|
19 |
virtual ~CSignalModel(){ delete model; }
|
20 |
RooAbsPdf *model;
|
21 |
};
|
22 |
|
23 |
class CBreitWignerConvCrystalBall : public CSignalModel
|
24 |
{
|
25 |
public:
|
26 |
CBreitWignerConvCrystalBall(RooRealVar &m, const Bool_t pass);
|
27 |
~CBreitWignerConvCrystalBall();
|
28 |
RooRealVar *mass, *width;
|
29 |
RooBreitWigner *bw;
|
30 |
RooRealVar *mean, *sigma, *alpha, *n;
|
31 |
RooCBShape *cb;
|
32 |
};
|
33 |
|
34 |
class CMCTemplateConvGaussian : public CSignalModel
|
35 |
{
|
36 |
public:
|
37 |
CMCTemplateConvGaussian(RooRealVar &m, TH1D* hist, const Bool_t pass, RooRealVar *sigma0=0, int intOrder=1);
|
38 |
~CMCTemplateConvGaussian();
|
39 |
RooRealVar *mean, *sigma;
|
40 |
RooGaussian *gaus;
|
41 |
TH1D *inHist;
|
42 |
RooDataHist *dataHist;
|
43 |
RooHistPdf *histPdf;
|
44 |
};
|
45 |
|
46 |
class CVoigtianCBShape : public CSignalModel
|
47 |
{
|
48 |
public:
|
49 |
CVoigtianCBShape(RooRealVar &m, const Bool_t pass);
|
50 |
~CVoigtianCBShape();
|
51 |
RooRealVar *mass, *width;
|
52 |
RooRealVar *sigma, *alpha, *n;
|
53 |
};
|
54 |
|
55 |
class CMCDatasetConvGaussian : public CSignalModel
|
56 |
{
|
57 |
public:
|
58 |
CMCDatasetConvGaussian(RooRealVar &m, TTree* tree, const Bool_t pass, RooRealVar *sigma0=0);
|
59 |
~CMCDatasetConvGaussian();
|
60 |
RooRealVar *mean, *sigma;
|
61 |
RooGaussian *gaus;
|
62 |
TTree *inTree;
|
63 |
RooDataSet *dataSet;
|
64 |
RooKeysPdf *keysPdf;
|
65 |
};
|
66 |
|
67 |
//--------------------------------------------------------------------------------------------------
|
68 |
CBreitWignerConvCrystalBall::CBreitWignerConvCrystalBall(RooRealVar &m, const Bool_t pass)
|
69 |
{
|
70 |
char name[10];
|
71 |
if(pass) sprintf(name,"%s","Pass");
|
72 |
else sprintf(name,"%s","Fail");
|
73 |
|
74 |
char vname[50];
|
75 |
|
76 |
sprintf(vname,"mass%s",name);
|
77 |
mass = new RooRealVar(vname,vname,91,80,100);
|
78 |
mass->setVal(91.1876);
|
79 |
mass->setConstant(kTRUE);
|
80 |
|
81 |
sprintf(vname,"width%s",name);
|
82 |
width = new RooRealVar(vname,vname,2.5,0.1,10);
|
83 |
width->setVal(2.4952);
|
84 |
width->setConstant(kTRUE);
|
85 |
|
86 |
sprintf(vname,"bw%s",name);
|
87 |
bw = new RooBreitWigner(vname,vname,m,*mass,*width);
|
88 |
|
89 |
if(pass) {
|
90 |
sprintf(vname,"mean%s",name); mean = new RooRealVar(vname,vname,0,-10,10);
|
91 |
sprintf(vname,"sigma%s",name); sigma = new RooRealVar(vname,vname,1,0.1,5);
|
92 |
sprintf(vname,"alpha%s",name); alpha = new RooRealVar(vname,vname,5,0,20);
|
93 |
sprintf(vname,"n%s",name); n = new RooRealVar(vname,vname,1,0,10);
|
94 |
} else {
|
95 |
sprintf(vname,"mean%s",name); mean = new RooRealVar(vname,vname,0,-10,10);
|
96 |
sprintf(vname,"sigma%s",name); sigma = new RooRealVar(vname,vname,1,0.1,5);
|
97 |
sprintf(vname,"alpha%s",name); alpha = new RooRealVar(vname,vname,5,0,20);
|
98 |
sprintf(vname,"n%s",name); n = new RooRealVar(vname,vname,1,0,10);
|
99 |
}
|
100 |
// n->setVal(1.0);
|
101 |
// n->setConstant(kTRUE);
|
102 |
|
103 |
sprintf(vname,"cb%s",name);
|
104 |
cb = new RooCBShape(vname,vname,m,*mean,*sigma,*alpha,*n);
|
105 |
|
106 |
sprintf(vname,"signal%s",name);
|
107 |
model = new RooFFTConvPdf(vname,vname,m,*bw,*cb);
|
108 |
}
|
109 |
|
110 |
CBreitWignerConvCrystalBall::~CBreitWignerConvCrystalBall()
|
111 |
{
|
112 |
delete mass;
|
113 |
delete width;
|
114 |
delete bw;
|
115 |
delete mean;
|
116 |
delete sigma;
|
117 |
delete alpha;
|
118 |
delete n;
|
119 |
delete cb;
|
120 |
}
|
121 |
|
122 |
//--------------------------------------------------------------------------------------------------
|
123 |
CMCTemplateConvGaussian::CMCTemplateConvGaussian(RooRealVar &m, TH1D* hist, const Bool_t pass, RooRealVar *sigma0, int intOrder)
|
124 |
{
|
125 |
char name[10];
|
126 |
if(pass) sprintf(name,"%s","Pass");
|
127 |
else sprintf(name,"%s","Fail");
|
128 |
|
129 |
char vname[50];
|
130 |
|
131 |
if(pass) {
|
132 |
sprintf(vname,"mean%s",name); mean = new RooRealVar(vname,vname,0,-10,10);
|
133 |
if(sigma0) { sigma = sigma0; }
|
134 |
else { sprintf(vname,"sigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
135 |
sprintf(vname,"gaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
136 |
} else {
|
137 |
sprintf(vname,"mean%s",name); mean = new RooRealVar(vname,vname,0,-10,10);
|
138 |
if(sigma0) { sigma = sigma0; }
|
139 |
else { sprintf(vname,"sigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
140 |
sprintf(vname,"gaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
141 |
}
|
142 |
|
143 |
|
144 |
sprintf(vname,"inHist_%s",hist->GetName());
|
145 |
inHist = (TH1D*)hist->Clone(vname);
|
146 |
|
147 |
sprintf(vname,"dataHist%s",name); dataHist = new RooDataHist(vname,vname,RooArgSet(m),inHist);
|
148 |
sprintf(vname,"histPdf%s",name); histPdf = new RooHistPdf(vname,vname,m,*dataHist,intOrder);
|
149 |
sprintf(vname,"signal%s",name); model = new RooFFTConvPdf(vname,vname,m,*histPdf,*gaus);
|
150 |
}
|
151 |
|
152 |
CMCTemplateConvGaussian::~CMCTemplateConvGaussian()
|
153 |
{
|
154 |
delete mean;
|
155 |
//delete sigma;
|
156 |
delete gaus;
|
157 |
delete inHist;
|
158 |
delete dataHist;
|
159 |
delete histPdf;
|
160 |
}
|
161 |
|
162 |
//--------------------------------------------------------------------------------------------------
|
163 |
CVoigtianCBShape::CVoigtianCBShape(RooRealVar &m, const Bool_t pass)
|
164 |
{
|
165 |
char name[10];
|
166 |
if(pass) sprintf(name,"%s","Pass");
|
167 |
else sprintf(name,"%s","Fail");
|
168 |
|
169 |
char vname[50];
|
170 |
|
171 |
sprintf(vname,"mass%s",name);
|
172 |
mass = new RooRealVar(vname,vname,91,80,100);
|
173 |
|
174 |
sprintf(vname,"width%s",name);
|
175 |
width = new RooRealVar(vname,vname,2.5,0.1,10);
|
176 |
width->setVal(2.4952);
|
177 |
width->setConstant(kTRUE);
|
178 |
|
179 |
if(pass) {
|
180 |
sprintf(vname,"sigma%s",name); sigma = new RooRealVar(vname,vname,1,0.1,10);
|
181 |
sprintf(vname,"alpha%s",name); alpha = new RooRealVar(vname,vname,5,0,20);
|
182 |
} else {
|
183 |
sprintf(vname,"sigma%s",name); sigma = new RooRealVar(vname,vname,1,0.1,10);
|
184 |
sprintf(vname,"alpha%s",name); alpha = new RooRealVar(vname,vname,5,0,20);
|
185 |
}
|
186 |
|
187 |
sprintf(vname,"n%s",name);
|
188 |
n = new RooRealVar(vname,vname,1,0,10);
|
189 |
n->setVal(1.0);
|
190 |
|
191 |
sprintf(vname,"signal%s",name);
|
192 |
model = new RooVoigtianShape(vname,vname,m,*mass,*sigma,*alpha,*n,*width,0);
|
193 |
}
|
194 |
|
195 |
CVoigtianCBShape::~CVoigtianCBShape()
|
196 |
{
|
197 |
delete mass;
|
198 |
delete width;
|
199 |
delete sigma;
|
200 |
delete alpha;
|
201 |
delete n;
|
202 |
}
|
203 |
|
204 |
//--------------------------------------------------------------------------------------------------
|
205 |
CMCDatasetConvGaussian::CMCDatasetConvGaussian(RooRealVar &m, TTree* tree, const Bool_t pass, RooRealVar *sigma0)
|
206 |
{
|
207 |
char name[10];
|
208 |
if(pass) sprintf(name,"%s","Pass");
|
209 |
else sprintf(name,"%s","Fail");
|
210 |
|
211 |
char vname[50];
|
212 |
|
213 |
if(pass) {
|
214 |
sprintf(vname,"mean%s",name); mean = new RooRealVar(vname,vname,0,-10,10);
|
215 |
if(sigma0) { sigma = sigma0; }
|
216 |
else { sprintf(vname,"sigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
217 |
sprintf(vname,"gaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
218 |
} else {
|
219 |
sprintf(vname,"mean%s",name); mean = new RooRealVar(vname,vname,0,-10,10);
|
220 |
if(sigma0) { sigma = sigma0; }
|
221 |
else { sprintf(vname,"sigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
222 |
sprintf(vname,"gaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
223 |
}
|
224 |
|
225 |
sprintf(vname,"inTree_%s",tree->GetName());
|
226 |
inTree = (TTree*)tree->Clone(vname);
|
227 |
|
228 |
sprintf(vname,"dataSet%s",name); dataSet = new RooDataSet(vname,vname,inTree,RooArgSet(m));
|
229 |
sprintf(vname,"keysPdf%s",name); keysPdf = new RooKeysPdf(vname,vname,m,*dataSet,RooKeysPdf::NoMirror,1);
|
230 |
sprintf(vname,"signal%s",name); model = new RooFFTConvPdf(vname,vname,m,*keysPdf,*gaus);
|
231 |
}
|
232 |
|
233 |
CMCDatasetConvGaussian::~CMCDatasetConvGaussian()
|
234 |
{
|
235 |
delete mean;
|
236 |
delete sigma;
|
237 |
delete gaus;
|
238 |
delete inTree;
|
239 |
delete dataSet;
|
240 |
delete keysPdf;
|
241 |
}
|