Logo ROOT   6.12/07
Reference Guide
RooBukinPdf.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * RW, Ruddick William UC Colorado wor@slac.stanford.edu *
7  * *
8  * Copyright (c) 2000-2005, Regents of the University of California *
9  * and Stanford University. All rights reserved. *
10  * *
11  * Redistribution and use in source and binary forms, *
12  * with or without modification, are permitted according to the terms *
13  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
14  *****************************************************************************/
15 
16 /** \class RooBukinPdf
17  \ingroup Roofit
18 
19 RooBukinPdf implements the NovosibirskA function
20 
21 #### Original Fortran Header below
22 ~~~
23  * Fitting function for asymmetric peaks with 6 free parameters:
24  * Ap - peak value
25  * Xp - peak position
26  * sigp - FWHM divided by 2*sqrt(2*log(2))=2.35
27  * xi - peak asymmetry parameter
28  * rho1 - parameter of the "left tail"
29  * rho2 - parameter of the "right tail"
30  * ---------------------------------------------
31  * May 26, 2003
32  * A.Bukin, Budker INP, Novosibirsk
33  * Documentation:
34  * http://www.slac.stanford.edu/BFROOT/www/Organization/CollabMtgs/2003/detJuly2003/Tues3a/bukin.ps
35  * -------------------------------------------
36 ~~~
37 **/
38 
39 #include "RooFit.h"
40 
41 #include <math.h>
42 
43 
44 #include "RooBukinPdf.h"
45 #include "RooRealVar.h"
46 #include "TMath.h"
47 
48 using namespace std;
49 
51 
52 ////////////////////////////////////////////////////////////////////////////////
53 
54 RooBukinPdf::RooBukinPdf(const char *name, const char *title,
55  RooAbsReal& _x, RooAbsReal& _Xp,
56  RooAbsReal& _sigp, RooAbsReal& _xi,
57  RooAbsReal& _rho1, RooAbsReal& _rho2) :
58  // The two addresses refer to our first dependent variable and
59  // parameter, respectively, as declared in the rdl file
60  RooAbsPdf(name, title),
61  x("x","x",this,_x),
62  Xp("Xp","Xp",this,_Xp),
63  sigp("sigp","sigp",this,_sigp),
64  xi("xi","xi",this,_xi),
65  rho1("rho1","rho1",this,_rho1),
66  rho2("rho2","rho2",this,_rho2)
67 {
68  // Constructor
69  consts = 2*sqrt(2*log(2.));
70 }
71 
72 ////////////////////////////////////////////////////////////////////////////////
73 
74 RooBukinPdf::RooBukinPdf(const RooBukinPdf& other, const char *name):
75  RooAbsPdf(other,name),
76  x("x",this,other.x),
77  Xp("Xp",this,other.Xp),
78  sigp("sigp",this,other.sigp),
79  xi("xi",this,other.xi),
80  rho1("rho1",this,other.rho1),
81  rho2("rho2",this,other.rho2)
82 
83 {
84  // Copy constructor
85  consts = 2*sqrt(2*log(2.));
86 }
87 
88 ////////////////////////////////////////////////////////////////////////////////
89 /// Implementation
90 
92 {
93  double r1=0,r2=0,r3=0,r4=0,r5=0,hp=0;
94  double x1 = 0,x2 = 0;
95  double fit_result = 0;
96 
97  hp=sigp*consts;
98  r3=log(2.);
99  r4=sqrt(TMath::Power(xi,2)+1);
100  r1=xi/r4;
101 
102  if(TMath::Abs(xi) > exp(-6.)){
103  r5=xi/log(r4+xi);
104  }
105  else
106  r5=1;
107 
108  x1 = Xp + (hp / 2) * (r1-1);
109  x2 = Xp + (hp / 2) * (r1+1);
110 
111  //--- Left Side
112  if(x < x1){
113  r2=rho1*TMath::Power((x-x1)/(Xp-x1),2)-r3 + 4 * r3 * (x-x1)/hp * r5 * r4/TMath::Power((r4-xi),2);
114  }
115 
116 
117  //--- Center
118  else if(x < x2) {
119  if(TMath::Abs(xi) > exp(-6.)) {
120  r2=log(1 + 4 * xi * r4 * (x-Xp)/hp)/log(1+2*xi*(xi-r4));
121  r2=-r3*(TMath::Power(r2,2));
122  }
123  else{
124  r2=-4*r3*TMath::Power(((x-Xp)/hp),2);
125  }
126  }
127 
128 
129  //--- Right Side
130  else {
131  r2=rho2*TMath::Power((x-x2)/(Xp-x2),2)-r3 - 4 * r3 * (x-x2)/hp * r5 * r4/TMath::Power((r4+xi),2);
132  }
133 
134  if(TMath::Abs(r2) > 100){
135  fit_result = 0;
136  }
137  else{
138  //---- Normalize the result
139  fit_result = exp(r2);
140  }
141 
142  return fit_result;
143 
144 }
RooRealProxy rho2
Definition: RooBukinPdf.h:49
RooRealProxy rho1
Definition: RooBukinPdf.h:48
RooRealProxy x
Definition: RooBukinPdf.h:44
STL namespace.
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:627
double sqrt(double)
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
RooRealProxy sigp
Definition: RooBukinPdf.h:46
RooRealProxy xi
Definition: RooBukinPdf.h:47
double consts
Definition: RooBukinPdf.h:55
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:359
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooRealProxy Xp
Definition: RooBukinPdf.h:45
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Double_t evaluate() const
Implementation.
Definition: RooBukinPdf.cxx:91
RooBukinPdf implements the NovosibirskA function.
Definition: RooBukinPdf.h:29
double exp(double)
char name[80]
Definition: TGX11.cxx:109
double log(double)