46 ClassImp(RooBernstein);
50 RooBernstein::RooBernstein()
57 RooBernstein::RooBernstein(
const char* name,
const char* title,
58 RooAbsReal& x,
const RooArgList& coefList):
59 RooAbsPdf(name, title),
60 _x(
"x",
"Dependent", this, x),
61 _coefList(
"coefficients",
"List of coefficients",this)
63 TIterator* coefIter = coefList.createIterator() ;
65 while((coef = (RooAbsArg*)coefIter->Next())) {
66 if (!dynamic_cast<RooAbsReal*>(coef)) {
67 cout <<
"RooBernstein::ctor(" << GetName() <<
") ERROR: coefficient " << coef->GetName()
68 <<
" is not of type RooAbsReal" << endl ;
71 _coefList.add(*coef) ;
78 RooBernstein::RooBernstein(
const RooBernstein& other,
const char* name) :
79 RooAbsPdf(other, name),
80 _x(
"x", this, other._x),
81 _coefList(
"coefList",this,other._coefList)
87 Double_t RooBernstein::evaluate()
const
89 Double_t xmin = _x.min();
90 Double_t x = (_x - xmin) / (_x.max() - xmin);
91 Int_t degree = _coefList.getSize() - 1;
92 RooFIter iter = _coefList.fwdIterator();
96 return ((RooAbsReal *)iter.next())->getVal();
98 }
else if(degree == 1) {
100 Double_t a0 = ((RooAbsReal *)iter.next())->getVal();
101 Double_t a1 = ((RooAbsReal *)iter.next())->getVal() - a0;
104 }
else if(degree == 2) {
106 Double_t a0 = ((RooAbsReal *)iter.next())->getVal();
107 Double_t a1 = 2 * (((RooAbsReal *)iter.next())->getVal() - a0);
108 Double_t a2 = ((RooAbsReal *)iter.next())->getVal() - a1 - a0;
109 return (a2 * x + a1) * x + a0;
111 }
else if(degree > 2) {
116 Double_t result = ((RooAbsReal *)iter.next())->getVal() * s;
117 for(Int_t i = 1; i < degree; i++) {
118 result = (result + t * TMath::Binomial(degree, i) * ((RooAbsReal *)iter.next())->getVal()) * s;
121 result += t * ((RooAbsReal *)iter.next())->getVal();
127 return TMath::SignalingNaN();
135 void compute(
size_t batchSize,
double xmax,
double xmin,
136 double * __restrict output,
137 const double * __restrict
const xData,
138 const RooListProxy& coefList)
140 constexpr
size_t block = 128;
141 const int nCoef = coefList.size();
142 const int degree = nCoef-1;
143 double X[block], _1_X[block], powX[block], pow_1_X[block];
144 double *Binomial =
new double[nCoef+5];
148 for (
int i=1; i<=degree; i++) {
149 Binomial[i] = Binomial[i-1]*(degree-i+1)/i;
152 for (
size_t i=0; i<batchSize; i+=block) {
153 const size_t stop = (i+block > batchSize) ? batchSize-i : block;
156 for (
size_t j=0; j<stop; j++) {
157 powX[j] = pow_1_X[j] = 1.0;
158 X[j] = (xData[i+j]-xmin) / (xmax-xmin);
164 for (
int k=2; k<=degree; k+=2)
165 for (
size_t j=0; j<stop; j++)
166 pow_1_X[j] *= _1_X[j]*_1_X[j];
169 for (
size_t j=0; j<stop; j++)
170 pow_1_X[j] *= _1_X[j];
173 for (
size_t j=0; j<stop; j++)
176 for (
int k=0; k<nCoef; k++) {
177 double coef =
static_cast<RooAbsReal&
>(coefList[k]).getVal();
178 for (
size_t j=0; j<stop; j++) {
179 output[i+j] += coef*Binomial[k]*powX[j]*pow_1_X[j];
183 pow_1_X[j] *= _1_X[j];
193 RooSpan<double> RooBernstein::evaluateBatch(std::size_t begin, std::size_t batchSize)
const {
194 auto xData = _x.getValBatch(begin, batchSize);
199 batchSize = xData.size();
200 auto output = _batchData.makeWritableBatchUnInit(begin, batchSize);
201 const double xmax = _x.max();
202 const double xmin = _x.min();
203 compute(batchSize, xmax, xmin, output.data(), xData.data(), _coefList);
210 Int_t RooBernstein::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars,
const char* rangeName)
const
212 if (rangeName && strlen(rangeName)) {
216 if (matchArgs(allVars, analVars, _x))
return 1;
222 Double_t RooBernstein::analyticalIntegral(Int_t code,
const char* rangeName)
const
225 Double_t xmin = _x.min(rangeName); Double_t xmax = _x.max(rangeName);
226 Int_t degree= _coefList.getSize()-1;
229 RooFIter iter = _coefList.fwdIterator() ;
231 for (
int i=0; i<=degree; ++i){
236 for (
int j=i; j<=degree; ++j){
237 temp += pow(-1.,j-i) * TMath::Binomial(degree, j) * TMath::Binomial(j,i) / (j+1);
239 temp *= ((RooAbsReal*)iter.next())->getVal();