00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031 #include "RooRarFit/rarVersion.hh"
00032
00033 #include "Riostream.h"
00034
00035 #include "RooRarFit/RooFlatte.hh"
00036 #include "RooFitCore/RooAbsReal.hh"
00037 #include "RooFitCore/RooRealVar.hh"
00038 #include "RooFitCore/RooComplex.hh"
00039
00040 ClassImp(RooFlatte)
00041
00042
00043 RooFlatte::RooFlatte(const char *name, const char *title,
00044 RooAbsReal& _x, RooAbsReal& _mean,
00045 RooAbsReal& _g0, RooAbsReal& _m0a, RooAbsReal& _m0b,
00046 RooAbsReal& _g1, RooAbsReal& _m1a, RooAbsReal& _m1b) :
00047 RooAbsPdf(name,title),
00048 x("x","Dependent",this,_x),
00049 mean("mean","Mean",this,_mean),
00050 g0("g0","Channel 1 coupling",this,_g0),
00051 m0a("m0a","Mass of particle 1 in channel 1",this,_m0a),
00052 m0b("m0b","Mass of particle 2 in channel 1",this,_m0b),
00053 g1("g1","Channel 2 coupling",this,_g1),
00054 m1a("m1a","Mass of particle 1 in channel 2",this,_m1a),
00055 m1b("m1b","Mass of particle 2 in channel 2",this,_m1b)
00056 {
00057 }
00058
00059
00060 RooFlatte::RooFlatte(const RooFlatte& other,
00061 const char* name) :
00062 RooAbsPdf(other,name),
00063 x("x",this,other.x),
00064 mean("mean",this,other.mean),
00065 g0("g0",this,other.g0),
00066 m0a("m0a",this,other.m0a),
00067 m0b("m0b",this,other.m0b),
00068 g1("g1",this,other.g1),
00069 m1a("m1a",this,other.m1a),
00070 m1b("m1b",this,other.m1b)
00071 {
00072 }
00073
00074
00075 Double_t RooFlatte::evaluate() const
00076 {
00077
00078
00079 if (g0<0 || g1<0) {return(0);}
00080
00081 Double_t s = x*x;
00082
00083
00084 Double_t E0a = 0.5 * (s + m0a*m0a - m0b*m0b) / x;
00085 Double_t qSq0 = E0a*E0a - m0a*m0a;
00086
00087
00088 Double_t E1a = 0.5 * (s + m1a*m1a - m1b*m1b) / x;
00089 Double_t qSq1 = E1a*E1a - m1a*m1a;
00090
00091 RooComplex gamma0 = (qSq0 > 0) ? RooComplex(g0*sqrt(qSq0),0) :
00092 RooComplex(0, g0*sqrt(-qSq0));
00093
00094 RooComplex gamma1 = (qSq1 > 0) ? RooComplex(g1*sqrt(qSq1),0) :
00095 RooComplex(0, g1*sqrt(-qSq1));
00096
00097 RooComplex gamma = gamma0 + gamma1;
00098
00099 RooComplex partB = RooComplex(0.0, 2*mean/x) * gamma;
00100 RooComplex partA(mean*mean - s, 0);
00101
00102 RooComplex denom = partA - partB;
00103
00104
00105 RooComplex T(1,0);
00106 T = T / denom;
00107
00108 return(T.abs2());
00109
00110 }
00111