40 ClassImp(RooBukinPdf);
 
   52 RooBukinPdf::RooBukinPdf(
const char *name, 
const char *title,
 
   53           RooAbsReal& _x,    RooAbsReal& _Xp,
 
   54           RooAbsReal& _sigp, RooAbsReal& _xi,
 
   55           RooAbsReal& _rho1, RooAbsReal& _rho2) :
 
   58   RooAbsPdf(name, title),
 
   60   Xp(
"Xp",
"Xp",this,_Xp),
 
   61   sigp(
"sigp",
"sigp",this,_sigp),
 
   62   xi(
"xi",
"xi",this,_xi),
 
   63   rho1(
"rho1",
"rho1",this,_rho1),
 
   64   rho2(
"rho2",
"rho2",this,_rho2)
 
   66     RooHelpers::checkRangeOfParameters(
this, {&_sigp}, 0.0);
 
   67     RooHelpers::checkRangeOfParameters(
this, {&_rho1},-1.0, 0.0);
 
   68     RooHelpers::checkRangeOfParameters(
this, {&_rho2}, 0.0, 1.0);
 
   69     RooHelpers::checkRangeOfParameters(
this, {&_xi}, -1.0, 1.0);
 
   74 RooBukinPdf::RooBukinPdf(
const RooBukinPdf& other, 
const char *name):
 
   75   RooAbsPdf(other,name),
 
   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)
 
   89 Double_t RooBukinPdf::evaluate()
 const 
   91   const double consts = 2*sqrt(2*log(2.0));
 
   92   double r1=0,r2=0,r3=0,r4=0,r5=0,hp=0;
 
   94   double fit_result = 0;
 
  101   if(fabs(xi) > exp(-6.)){
 
  107   x1 = Xp + (hp / 2) * (r1-1);
 
  108   x2 = Xp + (hp / 2) * (r1+1);
 
  112     r2=rho1*(x-x1)*(x-x1)/(Xp-x1)/(Xp-x1)-r3 + 4 * r3 * (x-x1)/hp * r5 * r4/(r4-xi)/(r4-xi);
 
  118     if(fabs(xi) > exp(-6.)) {
 
  119       r2=log(1 + 4 * xi * r4 * (x-Xp)/hp)/log(1+2*xi*(xi-r4));
 
  123       r2=-4*r3*(x-Xp)*(x-Xp)/hp/hp;
 
  130     r2=rho2*(x-x2)*(x-x2)/(Xp-x2)/(Xp-x2)-r3 - 4 * r3 * (x-x2)/hp * r5 * r4/(r4+xi)/(r4+xi);
 
  138     fit_result = exp(r2);
 
  149 template<
class Tx, 
class TXp, 
class TSigp, 
class Txi, 
class Trho1, 
class Trho2>
 
  150 void compute(  
size_t batchSize,
 
  151                double * __restrict output,
 
  152                Tx X, TXp XP, TSigp SP, Txi XI, Trho1 R1, Trho2 R2)
 
  154   const double r3 = log(2.0);
 
  155   const double r6 = exp(-6.0);
 
  156   const double r7 = 2*sqrt(2*log(2.0));
 
  158   for (
size_t i=0; i<batchSize; i++) {
 
  159     const double r1 = XI[i]/sqrt(XI[i]*XI[i]+1);
 
  160     const double r4 = sqrt(XI[i]*XI[i]+1);
 
  161     const double hp = 1 / (SP[i]*r7);
 
  162     const double x1 = XP[i] + 0.5*SP[i]*r7*(r1-1);
 
  163     const double x2 = XP[i] + 0.5*SP[i]*r7*(r1+1);
 
  166     if (XI[i]>r6 || XI[i]<-r6) r5 = XI[i]/log(r4+XI[i]);
 
  168     double factor=1, y=X[i]-x1, Yp=XP[i]-x1, yi=r4-XI[i], rho=R1[i];
 
  177     output[i] = rho*y*y/Yp/Yp -r3 + factor*4*r3*y*hp*r5*r4/yi/yi;
 
  178     if (X[i]>=x1 && X[i]<x2) {
 
  179       output[i] = _rf_fast_log(1 + 4*XI[i]*r4*(X[i]-XP[i])*hp) / _rf_fast_log(1 +2*XI[i]*( XI[i]-r4 ));
 
  180       output[i] *= -output[i]*r3;
 
  182     if (X[i]>=x1 && X[i]<x2 && XI[i]<r6 && XI[i]>-r6) {
 
  183       output[i] = -4*r3*(X[i]-XP[i])*(X[i]-XP[i])*hp*hp;
 
  186   for (
size_t i=0; i<batchSize; i++) {
 
  187     output[i] = _rf_fast_exp(output[i]);
 
  193 RooSpan<double> RooBukinPdf::evaluateBatch(std::size_t begin, std::size_t batchSize)
 const {
 
  194   using namespace BatchHelpers;
 
  196   EvaluateInfo info = getInfo( {&x, &Xp, &sigp, &xi, &rho1, &rho2}, begin, batchSize );
 
  197   if (info.nBatches == 0) {
 
  200   auto output = _batchData.makeWritableBatchUnInit(begin, batchSize);
 
  201   auto xData = x.getValBatch(begin, info.size);
 
  203   if (info.nBatches==1 && !xData.empty()) {
 
  204     compute(info.size, output.data(), xData.data(),
 
  205     BracketAdapter<double> (Xp),
 
  206     BracketAdapter<double> (sigp),
 
  207     BracketAdapter<double> (xi),
 
  208     BracketAdapter<double> (rho1),
 
  209     BracketAdapter<double> (rho2));
 
  212     compute(info.size, output.data(),
 
  213     BracketAdapterWithMask (x,x.getValBatch(begin,info.size)),
 
  214     BracketAdapterWithMask (Xp,Xp.getValBatch(begin,info.size)),
 
  215     BracketAdapterWithMask (sigp,sigp.getValBatch(begin,info.size)),
 
  216     BracketAdapterWithMask (xi,xi.getValBatch(begin,info.size)),
 
  217     BracketAdapterWithMask (rho1,rho1.getValBatch(begin,info.size)),
 
  218     BracketAdapterWithMask (rho2,rho2.getValBatch(begin,info.size)));