1 |
#include "RooAbsPdf.h"
|
2 |
#include "RooAddPdf.h"
|
3 |
#include "RooRealVar.h"
|
4 |
#include "RooExponential.h"
|
5 |
#include "RooCMSShape.h"
|
6 |
#include "RooGenericPdf.h"
|
7 |
#include "RooChebychev.h"
|
8 |
#include "RooBernstein.h"
|
9 |
#include "RooPolynomial.h"
|
10 |
|
11 |
class CBackgroundModel
|
12 |
{
|
13 |
public:
|
14 |
CBackgroundModel():model(0){}
|
15 |
virtual ~CBackgroundModel() { delete model; }
|
16 |
RooAbsPdf *model;
|
17 |
};
|
18 |
|
19 |
class CExponential : public CBackgroundModel
|
20 |
{
|
21 |
public:
|
22 |
CExponential(RooRealVar &m, const Bool_t pass);
|
23 |
~CExponential();
|
24 |
RooRealVar *t;
|
25 |
};
|
26 |
|
27 |
class CErfExpo : public CBackgroundModel
|
28 |
{
|
29 |
public:
|
30 |
CErfExpo(RooRealVar &m, const Bool_t pass);
|
31 |
~CErfExpo();
|
32 |
RooRealVar *alfa, *beta, *gamma, *peak;
|
33 |
};
|
34 |
|
35 |
class CDoubleExp : public CBackgroundModel
|
36 |
{
|
37 |
public:
|
38 |
CDoubleExp(RooRealVar &m, const Bool_t pass);
|
39 |
~CDoubleExp();
|
40 |
RooExponential *exp1, *exp2;
|
41 |
RooRealVar *t1, *t2, *frac;
|
42 |
};
|
43 |
|
44 |
class CLinearExp : public CBackgroundModel
|
45 |
{
|
46 |
public:
|
47 |
CLinearExp(RooRealVar &m, const Bool_t pass);
|
48 |
~CLinearExp();
|
49 |
RooRealVar *a, *t;
|
50 |
};
|
51 |
|
52 |
class CQuadraticExp : public CBackgroundModel
|
53 |
{
|
54 |
public:
|
55 |
CQuadraticExp(RooRealVar &m, const Bool_t pass);
|
56 |
~CQuadraticExp();
|
57 |
RooRealVar *a1, *a2, *t;
|
58 |
};
|
59 |
|
60 |
class CQuadratic : public CBackgroundModel
|
61 |
{
|
62 |
public:
|
63 |
CQuadratic(RooRealVar &m, const Bool_t pass);
|
64 |
~CQuadratic();
|
65 |
RooRealVar *a1, *a2;
|
66 |
};
|
67 |
|
68 |
class CSecondOrderChebychev : public CBackgroundModel
|
69 |
{
|
70 |
public:
|
71 |
CSecondOrderChebychev(RooRealVar &m, const Bool_t pass);
|
72 |
~CSecondOrderChebychev();
|
73 |
RooRealVar *a1, *a2, *a3, *a4, *a5 ;
|
74 |
};
|
75 |
|
76 |
class CPolynomial : public CBackgroundModel
|
77 |
{
|
78 |
public:
|
79 |
CPolynomial(RooRealVar &m, const Bool_t pass);
|
80 |
~CPolynomial();
|
81 |
RooRealVar *a1, *a2, *a3, *a4, *a5;
|
82 |
};
|
83 |
|
84 |
|
85 |
class CMCBkgTemplate : public CBackgroundModel
|
86 |
{
|
87 |
public:
|
88 |
CMCBkgTemplate(RooRealVar &m, TH1D* hist, const Bool_t pass);
|
89 |
~CMCBkgTemplate();
|
90 |
TH1D *inHist;
|
91 |
RooDataHist *dataHist;
|
92 |
};
|
93 |
|
94 |
class CMCBkgTemplateConvGaussian : public CBackgroundModel
|
95 |
{
|
96 |
public:
|
97 |
CMCBkgTemplateConvGaussian(RooRealVar &m, TH1D* hist, const Bool_t pass,
|
98 |
RooRealVar *mean0 = 0, RooRealVar *sigma0=0);
|
99 |
~CMCBkgTemplateConvGaussian();
|
100 |
RooRealVar *mean, *sigma;
|
101 |
RooGaussian *gaus;
|
102 |
TH1D *inHist;
|
103 |
RooDataHist *dataHist;
|
104 |
RooHistPdf *histPdf;
|
105 |
};
|
106 |
|
107 |
|
108 |
class CMCBkgTemplateConvGaussianPlusExp : public CBackgroundModel
|
109 |
{
|
110 |
public:
|
111 |
CMCBkgTemplateConvGaussianPlusExp(RooRealVar &m, TH1D* hist, const Bool_t pass,
|
112 |
RooRealVar *mean0 = 0, RooRealVar *sigma0=0, string label="");
|
113 |
CMCBkgTemplateConvGaussianPlusExp(RooRealVar &m, RooHistPdf* myHistPdf, const Bool_t pass,
|
114 |
RooRealVar *mean0, RooRealVar *sigma0, string label="");
|
115 |
~CMCBkgTemplateConvGaussianPlusExp();
|
116 |
RooRealVar *mean, *sigma;
|
117 |
RooRealVar *t;
|
118 |
RooRealVar *frac;
|
119 |
RooGaussian *gaus;
|
120 |
RooExponential *exp;
|
121 |
TH1D *inHist;
|
122 |
RooDataHist *dataHist;
|
123 |
RooHistPdf *histPdf;
|
124 |
RooFFTConvPdf *convPdf;
|
125 |
|
126 |
|
127 |
};
|
128 |
|
129 |
|
130 |
//--------------------------------------------------------------------------------------------------
|
131 |
CExponential::CExponential(RooRealVar &m, const Bool_t pass)
|
132 |
{
|
133 |
char name[10];
|
134 |
if(pass) sprintf(name,"%s","Pass");
|
135 |
else sprintf(name,"%s","Fail");
|
136 |
|
137 |
char vname[50];
|
138 |
|
139 |
sprintf(vname,"t%s",name);
|
140 |
if(pass)
|
141 |
t = new RooRealVar(vname,vname,-0.1,-1.,0.);
|
142 |
else
|
143 |
t = new RooRealVar(vname,vname,-0.1,-1.,0.);
|
144 |
|
145 |
sprintf(vname,"background%s",name);
|
146 |
model = new RooExponential(vname,vname,m,*t);
|
147 |
}
|
148 |
|
149 |
CExponential::~CExponential()
|
150 |
{
|
151 |
delete t;
|
152 |
}
|
153 |
|
154 |
//--------------------------------------------------------------------------------------------------
|
155 |
CErfExpo::CErfExpo(RooRealVar &m, const Bool_t pass)
|
156 |
{
|
157 |
char name[10];
|
158 |
if(pass) sprintf(name,"%s","Pass");
|
159 |
else sprintf(name,"%s","Fail");
|
160 |
|
161 |
char vname[50];
|
162 |
|
163 |
if(pass) {
|
164 |
sprintf(vname,"alfa%s",name); alfa = new RooRealVar(vname,vname,50,5,200);
|
165 |
sprintf(vname,"beta%s",name); beta = new RooRealVar(vname,vname,0.01,0,10);
|
166 |
sprintf(vname,"gamma%s",name); gamma = new RooRealVar(vname,vname,0.1,0,1);
|
167 |
} else {
|
168 |
sprintf(vname,"alfa%s",name); alfa = new RooRealVar(vname,vname,50,5,200);
|
169 |
sprintf(vname,"beta%s",name); beta = new RooRealVar(vname,vname,0.01,0,10);
|
170 |
sprintf(vname,"gamma%s",name); gamma = new RooRealVar(vname,vname,0.1,0,1);
|
171 |
}
|
172 |
|
173 |
sprintf(vname,"peak%s",name);
|
174 |
peak = new RooRealVar(vname,vname,91.1876,85,97);
|
175 |
peak->setVal(91.1876);
|
176 |
peak->setConstant(kTRUE);
|
177 |
|
178 |
sprintf(vname,"background%s",name);
|
179 |
model = new RooCMSShape(vname,vname,m,*alfa,*beta,*gamma,*peak);
|
180 |
}
|
181 |
|
182 |
CErfExpo::~CErfExpo()
|
183 |
{
|
184 |
delete alfa;
|
185 |
delete beta;
|
186 |
delete gamma;
|
187 |
delete peak;
|
188 |
}
|
189 |
|
190 |
//--------------------------------------------------------------------------------------------------
|
191 |
CDoubleExp::CDoubleExp(RooRealVar &m, const Bool_t pass)
|
192 |
{
|
193 |
char name[10];
|
194 |
if(pass) sprintf(name,"%s","Pass");
|
195 |
else sprintf(name,"%s","Fail");
|
196 |
|
197 |
char vname[50];
|
198 |
|
199 |
if(pass) {
|
200 |
sprintf(vname,"t1%s",name); t1 = new RooRealVar(vname,vname,-0.20,-1.,0.);
|
201 |
sprintf(vname,"t2%s",name); t2 = new RooRealVar(vname,vname,-0.05,-1.,0.);
|
202 |
sprintf(vname,"frac%s",name); frac = new RooRealVar(vname,vname, 0.50, 0.,1.);
|
203 |
} else {
|
204 |
sprintf(vname,"t1%s",name); t1 = new RooRealVar(vname,vname,-0.20,-1.,0.);
|
205 |
sprintf(vname,"t2%s",name); t2 = new RooRealVar(vname,vname,-0.05,-1.,0.);
|
206 |
sprintf(vname,"frac%s",name); frac = new RooRealVar(vname,vname, 0.50, 0.,1.);
|
207 |
}
|
208 |
|
209 |
sprintf(vname,"exp1%s",name);
|
210 |
exp1 = new RooExponential(vname,vname,m,*t1);
|
211 |
sprintf(vname,"exp2%s",name);
|
212 |
exp2 = new RooExponential(vname,vname,m,*t2);
|
213 |
sprintf(vname,"background%s",name);
|
214 |
model = new RooAddPdf(vname,vname,RooArgList(*exp1,*exp2),RooArgList(*frac));
|
215 |
}
|
216 |
|
217 |
CDoubleExp::~CDoubleExp()
|
218 |
{
|
219 |
delete exp1;
|
220 |
delete exp2;
|
221 |
delete t1;
|
222 |
delete t2;
|
223 |
delete frac;
|
224 |
}
|
225 |
|
226 |
//--------------------------------------------------------------------------------------------------
|
227 |
CLinearExp::CLinearExp(RooRealVar &m, const Bool_t pass)
|
228 |
{
|
229 |
char name[10];
|
230 |
if(pass) sprintf(name,"%s","Pass");
|
231 |
else sprintf(name,"%s","Fail");
|
232 |
|
233 |
char aname[50];
|
234 |
sprintf(aname,"a%s",name);
|
235 |
a = new RooRealVar(aname,aname,-0,-10.,10.);
|
236 |
//a->setConstant(kTRUE);
|
237 |
|
238 |
char tname[50];
|
239 |
sprintf(tname,"t%s",name);
|
240 |
t = new RooRealVar(tname,tname,-1e-6,-10.,0.);
|
241 |
//t->setConstant(kTRUE);
|
242 |
|
243 |
char formula[200];
|
244 |
sprintf(formula,"(1+%s*m)*exp(%s*m)",aname,tname);
|
245 |
|
246 |
char vname[50]; sprintf(vname,"background%s",name);
|
247 |
model = new RooGenericPdf(vname,vname,formula,RooArgList(m,*a,*t));
|
248 |
}
|
249 |
|
250 |
CLinearExp::~CLinearExp()
|
251 |
{
|
252 |
delete a;
|
253 |
delete t;
|
254 |
}
|
255 |
|
256 |
//--------------------------------------------------------------------------------------------------
|
257 |
CQuadraticExp::CQuadraticExp(RooRealVar &m, const Bool_t pass)
|
258 |
{
|
259 |
char name[10];
|
260 |
if(pass) sprintf(name,"%s","Pass");
|
261 |
else sprintf(name,"%s","Fail");
|
262 |
|
263 |
char a1name[50];
|
264 |
sprintf(a1name,"a1%s",name);
|
265 |
a1 = new RooRealVar(a1name,a1name,0,-10,10.);
|
266 |
//a1->setConstant(kTRUE);
|
267 |
|
268 |
char a2name[50];
|
269 |
sprintf(a2name,"a2%s",name);
|
270 |
a2 = new RooRealVar(a2name,a2name,0.0,-10,10);
|
271 |
//a2->setConstant(kTRUE);
|
272 |
|
273 |
char tname[50];
|
274 |
sprintf(tname,"t%s",name);
|
275 |
t = new RooRealVar(tname,tname,-1e-6,-10.,0.);
|
276 |
//t->setConstant(kTRUE);
|
277 |
|
278 |
char formula[200];
|
279 |
sprintf(formula,"(1+%s*m+%s*m*m)*exp(%s*m)",a1name,a2name,tname);
|
280 |
|
281 |
char vname[50]; sprintf(vname,"background%s",name);
|
282 |
model = new RooGenericPdf(vname,vname,formula,RooArgList(m,*a1,*a2,*t));
|
283 |
}
|
284 |
|
285 |
CQuadraticExp::~CQuadraticExp()
|
286 |
{
|
287 |
delete a1;
|
288 |
delete a2;
|
289 |
delete t;
|
290 |
}
|
291 |
|
292 |
|
293 |
//--------------------------------------------------------------------------------------------------
|
294 |
CQuadratic::CQuadratic(RooRealVar &m, const Bool_t pass)
|
295 |
{
|
296 |
char name[10];
|
297 |
if(pass) sprintf(name,"%s","Pass");
|
298 |
else sprintf(name,"%s","Fail");
|
299 |
|
300 |
char a1name[50];
|
301 |
sprintf(a1name,"a1%s",name);
|
302 |
a1 = new RooRealVar(a1name,a1name,0,-10,10.);
|
303 |
//a1->setConstant(kTRUE);
|
304 |
|
305 |
char a2name[50];
|
306 |
sprintf(a2name,"a2%s",name);
|
307 |
a2 = new RooRealVar(a2name,a2name,0.0,-10,10);
|
308 |
//a2->setConstant(kTRUE);
|
309 |
|
310 |
char formula[200];
|
311 |
sprintf(formula,"(1+%s*m+%s*m*m)",a1name,a2name);
|
312 |
|
313 |
char vname[50]; sprintf(vname,"background%s",name);
|
314 |
model = new RooGenericPdf(vname,vname,formula,RooArgList(m,*a1,*a2));
|
315 |
}
|
316 |
|
317 |
CQuadratic::~CQuadratic()
|
318 |
{
|
319 |
delete a1;
|
320 |
delete a2;
|
321 |
}
|
322 |
|
323 |
|
324 |
|
325 |
//--------------------------------------------------------------------------------------------------
|
326 |
CSecondOrderChebychev::CSecondOrderChebychev(RooRealVar &m, const Bool_t pass)
|
327 |
{
|
328 |
char name[10];
|
329 |
if(pass) sprintf(name,"%s","Pass");
|
330 |
else sprintf(name,"%s","Fail");
|
331 |
|
332 |
char a1name[50];
|
333 |
sprintf(a1name,"a1%s",name);
|
334 |
a1 = new RooRealVar(a1name,a1name,0,-10,10.);
|
335 |
//a1->setConstant(kTRUE);
|
336 |
|
337 |
char a2name[50];
|
338 |
sprintf(a2name,"a2%s",name);
|
339 |
a2 = new RooRealVar(a2name,a2name,0.0,-10,10);
|
340 |
//a2->setConstant(kTRUE);
|
341 |
|
342 |
char a3name[50];
|
343 |
sprintf(a3name,"a3%s",name);
|
344 |
a3 = new RooRealVar(a3name,a3name,0.0,-10,10);
|
345 |
//a2->setConstant(kTRUE);
|
346 |
|
347 |
char a4name[50];
|
348 |
sprintf(a4name,"a4%s",name);
|
349 |
a4 = new RooRealVar(a4name,a4name,0.0,-10,10);
|
350 |
//a4->setConstant(kTRUE);
|
351 |
|
352 |
char a5name[50];
|
353 |
sprintf(a5name,"a5%s",name);
|
354 |
a5 = new RooRealVar(a5name,a5name,0.0,-10,10);
|
355 |
//a5->setConstant(kTRUE);
|
356 |
|
357 |
char vname[50]; sprintf(vname,"background%s",name);
|
358 |
model = new RooChebychev(vname,vname,m, RooArgSet(*a1,*a2, *a3, *a4, *a5));
|
359 |
|
360 |
}
|
361 |
|
362 |
CSecondOrderChebychev::~CSecondOrderChebychev()
|
363 |
{
|
364 |
delete a1;
|
365 |
delete a2;
|
366 |
}
|
367 |
|
368 |
|
369 |
//--------------------------------------------------------------------------------------------------
|
370 |
CPolynomial::CPolynomial(RooRealVar &m, const Bool_t pass)
|
371 |
{
|
372 |
char name[10];
|
373 |
if(pass) sprintf(name,"%s","Pass");
|
374 |
else sprintf(name,"%s","Fail");
|
375 |
|
376 |
char a1name[50];
|
377 |
sprintf(a1name,"a1%s",name);
|
378 |
a1 = new RooRealVar(a1name,a1name,0.027,-10, 10);
|
379 |
a1->setConstant(kTRUE);
|
380 |
|
381 |
char a2name[50];
|
382 |
sprintf(a2name,"a2%s",name);
|
383 |
a2 = new RooRealVar(a2name,a2name,-1.67e-4,-10, 10);
|
384 |
a2->setConstant(kTRUE);
|
385 |
|
386 |
char a3name[50];
|
387 |
sprintf(a3name,"a3%s",name);
|
388 |
a3 = new RooRealVar(a3name,a3name,0.0,-10,10);
|
389 |
//a3->setConstant(kTRUE);
|
390 |
|
391 |
char a4name[50];
|
392 |
sprintf(a4name,"a4%s",name);
|
393 |
a4 = new RooRealVar(a4name,a4name,0.0,-10,10);
|
394 |
//a4->setConstant(kTRUE);
|
395 |
|
396 |
char a5name[50];
|
397 |
sprintf(a5name,"a5%s",name);
|
398 |
a5 = new RooRealVar(a5name,a5name,0.0,-10,10);
|
399 |
//a5->setConstant(kTRUE);
|
400 |
|
401 |
char vname[50]; sprintf(vname,"background%s",name);
|
402 |
model = new RooPolynomial(vname,vname,m, RooArgSet(*a1,*a2));
|
403 |
|
404 |
}
|
405 |
|
406 |
CPolynomial::~CPolynomial()
|
407 |
{
|
408 |
delete a1;
|
409 |
delete a2;
|
410 |
// delete a3;
|
411 |
}
|
412 |
|
413 |
//--------------------------------------------------------------------------------------------------
|
414 |
CMCBkgTemplate::CMCBkgTemplate(RooRealVar &m, TH1D* hist, const Bool_t pass)
|
415 |
{
|
416 |
char name[10];
|
417 |
if(pass) sprintf(name,"%s","Pass");
|
418 |
else sprintf(name,"%s","Fail");
|
419 |
|
420 |
char vname[50];
|
421 |
|
422 |
sprintf(vname,"bkgInputHist_%s",hist->GetName());
|
423 |
inHist = (TH1D*)hist->Clone(vname);
|
424 |
|
425 |
sprintf(vname,"bkgDataHist%s",name); dataHist = new RooDataHist(vname,vname,RooArgSet(m),inHist);
|
426 |
sprintf(vname,"background%s",name); model = new RooHistPdf(vname,vname,m,*dataHist,8);
|
427 |
}
|
428 |
|
429 |
CMCBkgTemplate::~CMCBkgTemplate()
|
430 |
{
|
431 |
delete inHist;
|
432 |
delete dataHist;
|
433 |
}
|
434 |
|
435 |
//--------------------------------------------------------------------------------------------------
|
436 |
CMCBkgTemplateConvGaussian::CMCBkgTemplateConvGaussian(RooRealVar &m, TH1D* hist, const Bool_t pass,
|
437 |
RooRealVar *mean0, RooRealVar *sigma0)
|
438 |
{
|
439 |
char name[10];
|
440 |
if(pass) sprintf(name,"%s","Pass");
|
441 |
else sprintf(name,"%s","Fail");
|
442 |
|
443 |
char vname[50];
|
444 |
|
445 |
if(pass) {
|
446 |
if(mean0) { mean = mean0; }
|
447 |
else { sprintf(vname,"bkgmean%s",name); mean = new RooRealVar(vname,vname,0,-10,10); }
|
448 |
if(sigma0) { sigma = sigma0; }
|
449 |
else { sprintf(vname,"bkgsigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
450 |
sprintf(vname,"bkggaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
451 |
} else {
|
452 |
if(mean0) { mean = mean0; }
|
453 |
else { sprintf(vname,"bkgmean%s",name); mean = new RooRealVar(vname,vname,0,-10,10); }
|
454 |
if(sigma0) { sigma = sigma0; }
|
455 |
else { sprintf(vname,"bkgsigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
456 |
sprintf(vname,"bkggaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
457 |
}
|
458 |
|
459 |
sprintf(vname,"bkgInputHist_%s",hist->GetName());
|
460 |
inHist = (TH1D*)hist->Clone(vname);
|
461 |
|
462 |
sprintf(vname,"bkgDataHist%s",name); dataHist = new RooDataHist(vname,vname,RooArgSet(m),inHist);
|
463 |
sprintf(vname,"bkgHistPdf%s",name); histPdf = new RooHistPdf(vname,vname,m,*dataHist,8);
|
464 |
sprintf(vname,"background%s",name); model = new RooFFTConvPdf(vname,vname,m,*histPdf,*gaus);
|
465 |
}
|
466 |
|
467 |
CMCBkgTemplateConvGaussian::~CMCBkgTemplateConvGaussian()
|
468 |
{
|
469 |
delete inHist;
|
470 |
delete dataHist;
|
471 |
}
|
472 |
|
473 |
|
474 |
//--------------------------------------------------------------------------------------------------
|
475 |
CMCBkgTemplateConvGaussianPlusExp::CMCBkgTemplateConvGaussianPlusExp(RooRealVar &m, TH1D* hist, const Bool_t pass,
|
476 |
RooRealVar *mean0, RooRealVar *sigma0, string label)
|
477 |
{
|
478 |
char name[10];
|
479 |
if(pass) sprintf(name,"%s","Pass");
|
480 |
else sprintf(name,"%s","Fail");
|
481 |
|
482 |
char vname[50];
|
483 |
|
484 |
if(pass) {
|
485 |
if(mean0) { mean = mean0; }
|
486 |
else { sprintf(vname,"bkgmean%s",name); mean = new RooRealVar(vname,vname,0,-10,10); }
|
487 |
if(sigma0) { sigma = sigma0; }
|
488 |
else { sprintf(vname,"bkgsigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
489 |
sprintf(vname,"t%s",name); t = new RooRealVar(vname,vname,-0.20,-1.,0.);
|
490 |
sprintf(vname,"frac%s",name); frac = new RooRealVar(vname,vname, 0.20, 0.,1.);
|
491 |
sprintf(vname,"bkggaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
492 |
sprintf(vname,"bkgexp%s",name); exp = new RooExponential(vname,vname,m,*t);
|
493 |
} else {
|
494 |
if(mean0) { mean = mean0; }
|
495 |
else { sprintf(vname,"bkgmean%s",name); mean = new RooRealVar(vname,vname,0,-10,10); }
|
496 |
if(sigma0) { sigma = sigma0; }
|
497 |
else { sprintf(vname,"bkgsigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
498 |
sprintf(vname,"t%s",name); t = new RooRealVar(vname,vname,-0.20,-1.,0.);
|
499 |
sprintf(vname,"frac%s",name); frac = new RooRealVar(vname,vname, 0.20, 0.,1.);
|
500 |
sprintf(vname,"bkggaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
501 |
sprintf(vname,"bkgexp%s",name); exp = new RooExponential(vname,vname,m,*t);
|
502 |
}
|
503 |
|
504 |
sprintf(vname,"bkgInputHist_%s",hist->GetName());
|
505 |
inHist = (TH1D*)hist->Clone(vname);
|
506 |
|
507 |
sprintf(vname,"bkgDataHist%s",name); dataHist = new RooDataHist(vname,vname,RooArgSet(m),inHist);
|
508 |
sprintf(vname,"bkgHistPdf%s",name); histPdf = new RooHistPdf(vname,vname,m,*dataHist,8);
|
509 |
sprintf(vname,"bkgConvPdf%s",name); convPdf = new RooFFTConvPdf(vname,vname,m,*histPdf,*gaus);
|
510 |
|
511 |
sprintf(vname,"background%s",name);
|
512 |
model = new RooAddPdf(vname,vname,RooArgList(*convPdf,*exp),RooArgList(*frac));
|
513 |
|
514 |
|
515 |
}
|
516 |
|
517 |
//--------------------------------------------------------------------------------------------------
|
518 |
CMCBkgTemplateConvGaussianPlusExp::CMCBkgTemplateConvGaussianPlusExp(RooRealVar &m, RooHistPdf* myHistPdf, const Bool_t pass,
|
519 |
RooRealVar *mean0, RooRealVar *sigma0, string label)
|
520 |
{
|
521 |
string Label = ""; if (label != "") Label = "_"+label;
|
522 |
char name[10];
|
523 |
if(pass) sprintf(name,"%s%s","Pass",Label.c_str());
|
524 |
else sprintf(name,"%s%s","Fail",Label.c_str());
|
525 |
|
526 |
char vname[50];
|
527 |
|
528 |
if(pass) {
|
529 |
if(mean0) { mean = mean0; }
|
530 |
else { sprintf(vname,"bkgmean%s",name); mean = new RooRealVar(vname,vname,0,-10,10); }
|
531 |
if(sigma0) { sigma = sigma0; }
|
532 |
else { sprintf(vname,"bkgsigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
533 |
sprintf(vname,"t%s",name); t = new RooRealVar(vname,vname,-0.20,-1.,0.);
|
534 |
sprintf(vname,"frac%s",name); frac = new RooRealVar(vname,vname, 0.20, 0.,1.);
|
535 |
sprintf(vname,"bkggaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
536 |
sprintf(vname,"bkgexp%s",name); exp = new RooExponential(vname,vname,m,*t);
|
537 |
} else {
|
538 |
if(mean0) { mean = mean0; }
|
539 |
else { sprintf(vname,"bkgmean%s",name); mean = new RooRealVar(vname,vname,0,-10,10); }
|
540 |
if(sigma0) { sigma = sigma0; }
|
541 |
else { sprintf(vname,"bkgsigma%s",name); sigma = new RooRealVar(vname,vname,2,0,5); }
|
542 |
sprintf(vname,"t%s",name); t = new RooRealVar(vname,vname,-0.20,-1.,0.);
|
543 |
sprintf(vname,"frac%s",name); frac = new RooRealVar(vname,vname, 0.20, 0.,1.);
|
544 |
sprintf(vname,"bkggaus%s",name); gaus = new RooGaussian(vname,vname,m,*mean,*sigma);
|
545 |
sprintf(vname,"bkgexp%s",name); exp = new RooExponential(vname,vname,m,*t);
|
546 |
}
|
547 |
|
548 |
histPdf = myHistPdf;
|
549 |
sprintf(vname,"bkgConvPdf%s",name); convPdf = new RooFFTConvPdf(vname,vname,m,*histPdf,*gaus);
|
550 |
|
551 |
sprintf(vname,"background%s",name);
|
552 |
model = new RooAddPdf(vname,vname,RooArgList(*convPdf,*exp),RooArgList(*frac));
|
553 |
|
554 |
|
555 |
}
|
556 |
|
557 |
CMCBkgTemplateConvGaussianPlusExp::~CMCBkgTemplateConvGaussianPlusExp()
|
558 |
{
|
559 |
delete inHist;
|
560 |
delete dataHist;
|
561 |
delete histPdf;
|
562 |
delete convPdf;
|
563 |
}
|