ROOT  6.06/09
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 //////////////////////////////////////////////////////////////////////////////
17 //
18 // BEGIN_HTML
19 // RooBukinPdf implements the NovosibirskA function
20 // END_HTML
21 //
22 
23 // Original Fortran Header below
24 /*****************************************************************************
25  * Fitting function for asymmetric peaks with 6 free parameters: *
26  * Ap - peak value *
27  * Xp - peak position *
28  * sigp - FWHM divided by 2*sqrt(2*log(2))=2.35 *
29  * xi - peak asymmetry parameter *
30  * rho1 - parameter of the "left tail" *
31  * rho2 - parameter of the "right tail" *
32  * --------------------------------------------- *
33  * May 26, 2003 *
34  * A.Bukin, Budker INP, Novosibirsk *
35  * Documentation: *
36  * http://www.slac.stanford.edu/BFROOT/www/Organization/CollabMtgs/2003/detJuly2003/Tues3a/bukin.ps
37  * ------------------------------------------- *
38  *****************************************************************************/
39 
40 #include "RooFit.h"
41 
42 #include <math.h>
43 
44 
45 #include "RooBukinPdf.h"
46 #include "RooRealVar.h"
47 #include "TMath.h"
48 
49 using namespace std;
50 
52 
53 
54 
55 ////////////////////////////////////////////////////////////////////////////////
56 
57 RooBukinPdf::RooBukinPdf(const char *name, const char *title,
58  RooAbsReal& _x, RooAbsReal& _Xp,
59  RooAbsReal& _sigp, RooAbsReal& _xi,
60  RooAbsReal& _rho1, RooAbsReal& _rho2) :
61  // The two addresses refer to our first dependent variable and
62  // parameter, respectively, as declared in the rdl file
63  RooAbsPdf(name, title),
64  x("x","x",this,_x),
65  Xp("Xp","Xp",this,_Xp),
66  sigp("sigp","sigp",this,_sigp),
67  xi("xi","xi",this,_xi),
68  rho1("rho1","rho1",this,_rho1),
69  rho2("rho2","rho2",this,_rho2)
70 {
71  // Constructor
72  consts = 2*sqrt(2*log(2.));
73 }
74 
75 
76 
77 ////////////////////////////////////////////////////////////////////////////////
78 
79 RooBukinPdf::RooBukinPdf(const RooBukinPdf& other, const char *name):
80  RooAbsPdf(other,name),
81  x("x",this,other.x),
82  Xp("Xp",this,other.Xp),
83  sigp("sigp",this,other.sigp),
84  xi("xi",this,other.xi),
85  rho1("rho1",this,other.rho1),
86  rho2("rho2",this,other.rho2)
87 
88 {
89  // Copy constructor
90  consts = 2*sqrt(2*log(2.));
91 }
92 
93 
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 /// Implementation
97 
99 {
100  double r1=0,r2=0,r3=0,r4=0,r5=0,hp=0;
101  double x1 = 0,x2 = 0;
102  double fit_result = 0;
103 
104  hp=sigp*consts;
105  r3=log(2.);
106  r4=sqrt(TMath::Power(xi,2)+1);
107  r1=xi/r4;
108 
109  if(TMath::Abs(xi) > exp(-6.)){
110  r5=xi/log(r4+xi);
111  }
112  else
113  r5=1;
114 
115  x1 = Xp + (hp / 2) * (r1-1);
116  x2 = Xp + (hp / 2) * (r1+1);
117 
118  //--- Left Side
119  if(x < x1){
120  r2=rho1*TMath::Power((x-x1)/(Xp-x1),2)-r3 + 4 * r3 * (x-x1)/hp * r5 * r4/TMath::Power((r4-xi),2);
121  }
122 
123 
124  //--- Center
125  else if(x < x2) {
126  if(TMath::Abs(xi) > exp(-6.)) {
127  r2=log(1 + 4 * xi * r4 * (x-Xp)/hp)/log(1+2*xi*(xi-r4));
128  r2=-r3*(TMath::Power(r2,2));
129  }
130  else{
131  r2=-4*r3*TMath::Power(((x-Xp)/hp),2);
132  }
133  }
134 
135 
136  //--- Right Side
137  else {
138  r2=rho2*TMath::Power((x-x2)/(Xp-x2),2)-r3 - 4 * r3 * (x-x2)/hp * r5 * r4/TMath::Power((r4+xi),2);
139  }
140 
141  if(TMath::Abs(r2) > 100){
142  fit_result = 0;
143  }
144  else{
145  //---- Normalize the result
146  fit_result = exp(r2);
147  }
148 
149  return fit_result;
150 
151 }
RooRealProxy rho2
Definition: RooBukinPdf.h:65
RooRealProxy rho1
Definition: RooBukinPdf.h:64
RooRealProxy x
Definition: RooBukinPdf.h:60
STL namespace.
ClassImp(RooBukinPdf) RooBukinPdf
Definition: RooBukinPdf.cxx:51
Double_t evaluate() const
Implementation.
Definition: RooBukinPdf.cxx:98
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
double sqrt(double)
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
unsigned int r3[N_CITIES]
Definition: simanTSP.cxx:323
RooRealProxy sigp
Definition: RooBukinPdf.h:62
RooRealProxy xi
Definition: RooBukinPdf.h:63
unsigned int r1[N_CITIES]
Definition: simanTSP.cxx:321
static const double x1[5]
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
#define name(a, b)
Definition: linkTestLib0.cpp:5
RooRealProxy Xp
Definition: RooBukinPdf.h:61
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
double exp(double)
unsigned int r2[N_CITIES]
Definition: simanTSP.cxx:322
double log(double)