CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HZZ2L2QRooPdfs.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <math.h>
3 #include "TMath.h"
4 
5 //#include "../interface/RooDoubleCB.h"
6 //#include "../interface/RooFermi.h"
7 //#include "../interface/RooRelBW.h"
8 #include "../interface/HZZ2L2QRooPdfs.h"
9 #include "RooRealVar.h"
10 #include "RooRealConstant.h"
11 
12 using namespace RooFit;
13 ClassImp(RooCB)
14 
15  RooCB::RooCB(){}
16 
17 RooCB::RooCB(const char *name, const char *title,
18  RooAbsReal& _x,
19  RooAbsReal& _mean,
20  RooAbsReal& _width,
21  RooAbsReal& _alpha,
22  RooAbsReal& _n,
23  RooAbsReal& _theta
24  ) :
25  RooAbsPdf(name,title),
26  x("x","x",this,_x),
27  mean("mean","mean",this,_mean),
28  width("width","width",this,_width),
29  alpha("alpha","alpha",this,_alpha),
30  n("n","n",this,_n),
31  theta("theta","theta",this,_theta)
32 {
33 }
34 
35 RooCB::RooCB(const RooCB& other, const char* name) :
36  RooAbsPdf(other,name),
37  x("x",this,other.x),
38  mean("mean",this,other.mean),
39  width("width",this,other.width),
40  alpha("alpha",this,other.alpha),
41  n("n",this,other.n),
42  theta("theta",this,other.theta)
43 {
44 }
45 
46 double RooCB::evaluate() const
47 {
48  double a = cos(theta)*alpha - sin(theta)*width;
49  double w = sin(theta)*alpha + cos(theta)*width;
50 
51  double t = (x-mean)/w;
52  if(a<0) t = -t;
53 
54  double absa = fabs((double)a);
55 
56  double A = TMath::Power(n/absa,n)*exp(-0.5*absa*absa);
57  double B = n/absa-absa;
58 
59  if(t >= -absa){
60  return exp(-0.5*t*t);
61  }else{
62  return A/TMath::Power(B-t,n);
63  }
64 }
65 
66 
67  ClassImp(RooDoubleCB)
68 
70 
71  RooDoubleCB::RooDoubleCB(const char *name, const char *title,
72  RooAbsReal& _x,
73  RooAbsReal& _mean,
74  RooAbsReal& _width,
75  RooAbsReal& _alpha1,
76  RooAbsReal& _n1,
77  RooAbsReal& _alpha2,
78  RooAbsReal& _n2
79  ) :
80  RooAbsPdf(name,title),
81  x("x","x",this,_x),
82  mean("mean","mean",this,_mean),
83  width("width","width",this,_width),
84  alpha1("alpha1","alpha1",this,_alpha1),
85  n1("n1","n1",this,_n1),
86  alpha2("alpha2","alpha2",this,_alpha2),
87  n2("n2","n2",this,_n2)
88  {
89  }
90 
91 
92  RooDoubleCB::RooDoubleCB(const RooDoubleCB& other, const char* name) :
93  RooAbsPdf(other,name),
94  x("x",this,other.x),
95  mean("mean",this,other.mean),
96  width("width",this,other.width),
97  alpha1("alpha1",this,other.alpha1),
98  n1("n1",this,other.n1),
99  alpha2("alpha2",this,other.alpha2),
100  n2("n2",this,other.n2)
101 
102  {
103  }
104 
105  double RooDoubleCB::evaluate() const
106  {
107  double A1 = pow(n1/fabs(alpha1),n1)*exp(-alpha1*alpha1/2);
108  double A2 = pow(n2/fabs(alpha2),n2)*exp(-alpha2*alpha2/2);
109  double B1 = n1/fabs(alpha1)-fabs(alpha1);
110  double B2 = n2/fabs(alpha2)-fabs(alpha2);
111 
112  if((x-mean)/width>-alpha1 && (x-mean)/width<alpha2){
113  return exp(-(x-mean)*(x-mean)/(2*width*width));
114  }else if((x-mean)/width<-alpha1){
115  return A1*pow(B1-(x-mean)/width,-n1);
116  }else if((x-mean)/width>alpha2){
117  return A2*pow(B2+(x-mean)/width,-n2);
118  }else{
119  cout << "ERROR evaluating range..." << endl;
120  return 99;
121  }
122 
123  }
124 
125  ClassImp(RooFermi)
126 
128 
129  RooFermi::RooFermi(const char *name, const char *title,
130  RooAbsReal& _x,
131  RooAbsReal& _cutOff,
132  RooAbsReal& _beta
133  ) :
134  RooAbsPdf(name,title),
135  x("x","x",this,_x),
136  cutOff("cutOff","cutOff",this,_cutOff),
137  beta("beta","beta",this,_beta)
138  {
139  }
140 
141 
142  RooFermi::RooFermi(const RooFermi& other, const char* name) :
143  RooAbsPdf(other,name),
144  x("x",this,other.x),
145  cutOff("cutOff",this,other.cutOff),
146  beta("beta",this,other.beta)
147 
148  {
149  }
150 
151 
152 
153  double RooFermi::evaluate() const
154  {
155  return 1.0/(exp((cutOff-x)/beta)+1);
156  }
157 
158  ClassImp(RooRelBW)
159 
161 
162  RooRelBW::RooRelBW(const char *name, const char *title,
163  RooAbsReal& _x,
164  RooAbsReal& _mean,
165  RooAbsReal& _width,
166  RooAbsReal& _n
167  ) :
168  RooAbsPdf(name,title),
169  x("x","x",this,_x),
170  mean("mean","mean",this,_mean),
171  width("width","width",this,_width),
172  n("n","n",this,_n)
173  {
174  }
175 
176 
177  RooRelBW::RooRelBW(const RooRelBW& other, const char* name) :
178  RooAbsPdf(other,name),
179  x("x",this,other.x),
180  mean("mean",this,other.mean),
181  width("width",this,other.width),
182  n("n",this,other.n)
183 
184  {
185  }
186 
187 
188 
189  double RooRelBW::evaluate() const
190  {
191  return pow(x*x,n)/((x*x-mean*mean)*(x*x-mean*mean)+pow(x*x/(mean*mean),2*n)*mean*mean*width*width);
192  }
193 
194 
195 ClassImp(Triangle)
196 
198 
199 Triangle::Triangle(const char *name, const char *title,
200  RooAbsReal& _m,
201  RooAbsReal& _start,
202  RooAbsReal& _turn,
203  RooAbsReal& _stop
204  ):
205  RooAbsPdf(name, title),
206  m("m", "Dependent", this, _m),
207  start("start","start",this,_start),
208  turn("turn","turn",this,_turn),
209  stop("stop","stop",this,_stop)
210 {
211 }
212 
213 Triangle::Triangle(const Triangle& other, const char* name) :
214  RooAbsPdf(other, name), m("m", this, other.m),start("start", this, other.start), turn("turn", this, other.turn), stop("stop", this, other.stop)
215 {
216 }
217 
218 Double_t Triangle::evaluate() const
219 {
220  //std::cout << m << " "<<1.+(start-m)/turn << " " << 1+(turn-m)/stop << std::endl;
221  if(m<turn && m > turn+start)
222  return 1.+(turn-m)/start;
223  if(m>=turn && m < turn+stop)
224  return 1.+(turn-m)/stop;
225 
226  return 0;
227 }
228 
229 
230 Int_t Triangle::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* range) const
231 {
232  if (matchArgs(allVars,analVars,m)) return 1;
233  return 0;
234 }
235 
236 Double_t Triangle::analyticalIntegral(Int_t code, const char* rangeName) const
237 {
238 
239  // WARNING, ASSSUMES TURN TO BE IN INTERVAL
240  assert(code==1) ;
241  //whole triangle
242  Double_t sumleft = sqrt(1+ (turn+start)*(turn+start) ) ;
243  Double_t sumright= sqrt(1+ (turn+stop)*(turn+stop) );
244 
245 
246  if(m.min() < turn+start)// correct for left missing bit
247  sumleft -= sumleft*(m.min()-(turn+start))/fabs(start);
248 
249 
250  if(m.max() > turn+stop)// correct for right missing bit
251  sumright -= sumright*(turn+stop -m.max())/fabs(stop);
252 
253 
254 
255  return sumleft+sumright;
256 }
257 
258 
259 
260 
261 ClassImp(RooLevelledExp)
262 
264 
265 RooLevelledExp::RooLevelledExp(const char *name, const char *title,
266  RooAbsReal& _x,
267  RooAbsReal& _sigma,
268  RooAbsReal& _alpha,
269  RooAbsReal& _m,
270  RooAbsReal& _theta):
271  RooAbsPdf(name,title),
272  x("x","x",this,_x),
273  sigma("sigma","sigma",this,_sigma),
274  alpha("alpha","alpha",this,_alpha),
275  m("m","m",this,_m),
276  // k("k","k",this,_k),
277  theta("theta","theta",this,_theta)
278 {
279 }
280 
282  RooAbsPdf(other,name),
283  x("x",this,other.x),
284  sigma("sigma",this,other.sigma),
285  alpha("alpha",this,other.alpha),
286  m("m",this,other.m),
287  theta("theta",this,other.theta)
288 {
289 }
290 
292 {
293  double res=0.0;
294  double s = cos(theta)*sigma - sin(theta)*alpha;
295  double a = sin(theta)*sigma + cos(theta)*alpha;
296 
297  //original
298  double t = fabs(x-m);
299  double den = (s + a*t);
300  res=exp(-1.0*t/den);
301 
302 
303  return res;
304 }
RooRealProxy alpha2
const double beta
RooRealProxy m
float alpha
Definition: AMPTWrapper.h:95
RooRealProxy width
RooRealProxy sigma
RooRealProxy alpha1
RooRealProxy width
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
Double_t evaluate() const
RooRealProxy start
RooRealProxy m
Double_t evaluate() const
RooRealProxy alpha
RooRealProxy x
RooRealProxy n
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
RooRealProxy cutOff
Double_t evaluate() const
RooRealProxy beta
RooRealProxy theta
RooRealProxy x
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
RooRealProxy mean
Double_t evaluate() const
Double_t evaluate() const
RooRealProxy n
RooRealProxy n1
RooRealProxy x
RooRealProxy x
Double_t evaluate() const
Double_t mean(Double_t mH, TString proc)
double a
Definition: hdecay.h:121
RooRealProxy mean
RooRealProxy alpha
RooRealProxy stop
tuple cout
Definition: gather_cfg.py:121
RooRealProxy mean
RooRealProxy theta
x
Definition: VDTMath.h:216
RooRealProxy width
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
RooRealProxy turn
RooRealProxy x
T w() const
RooRealProxy n2