58 ClassImp(RooNonCentralChiSquare);
62 RooNonCentralChiSquare::RooNonCentralChiSquare(
const char *name,
const char *title,
65 RooAbsReal& _lambda) :
66 RooAbsPdf(name,title),
69 lambda(
"lambda",
"lambda",this,_lambda),
72 fHasIssuedConvWarning(false),
73 fHasIssuedSumWarning(false)
75 #ifdef R__HAS_MATHMORE
76 ccoutD(InputArguments) <<
"RooNonCentralChiSquare::ctor(" << GetName() <<
77 "MathMore Available, will use Bessel function expressions unless SetForceSum(true) "<< endl ;
86 RooNonCentralChiSquare::RooNonCentralChiSquare(
const RooNonCentralChiSquare& other,
const char* name) :
87 RooAbsPdf(other,name),
90 lambda(
"lambda",this,other.lambda),
91 fErrorTol(other.fErrorTol),
92 fMaxIters(other.fMaxIters),
93 fHasIssuedConvWarning(false),
94 fHasIssuedSumWarning(false)
96 #ifdef R__HAS_MATHMORE
97 ccoutD(InputArguments) <<
"RooNonCentralChiSquare::ctor(" << GetName() <<
98 "MathMore Available, will use Bessel function expressions unless SetForceSum(true) "<< endl ;
99 fForceSum = other.fForceSum;
107 void RooNonCentralChiSquare::SetForceSum(Bool_t flag) {
109 #ifndef R__HAS_MATHMORE
111 ccoutD(InputArguments) <<
"RooNonCentralChiSquare::SetForceSum" << GetName() <<
112 "MathMore is not available- ForceSum must be on "<< endl ;
121 Double_t RooNonCentralChiSquare::evaluate()
const
128 Double_t xmin = x.min();
129 Double_t xmax = x.max();
135 _x = xmin + 1e-3*(xmax-xmin);
140 return ROOT::Math::chisquared_pdf(_x,k);
149 if(!fHasIssuedSumWarning){
150 coutI(InputArguments) <<
"RooNonCentralChiSquare sum being forced" << endl ;
151 fHasIssuedSumWarning=
true;
155 double errorTol = fErrorTol;
156 int MaxIters = fMaxIters;
157 int iDominant = (int) TMath::Floor(lambda/2);
162 for(
int i = iDominant; ; ++i){
163 ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)*ROOT::Math::chisquared_pdf(_x,k+2*i)/TMath::Gamma(i+1);
166 if(ithTerm/sum < errorTol)
169 if( i>iDominant+MaxIters){
170 if(!fHasIssuedConvWarning){
171 fHasIssuedConvWarning=
true;
172 coutW(Eval) <<
"RooNonCentralChiSquare did not converge: for x=" << x <<
" k="<<k
173 <<
", lambda="<<lambda <<
" fractional error = " << ithTerm/sum
174 <<
"\n either adjust tolerance with SetErrorTolerance(tol) or max_iter with SetMaxIter(max_it)"
181 for(
int i = iDominant - 1; i >= 0; --i){
183 ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)*ROOT::Math::chisquared_pdf(_x,k+2*i)/TMath::Gamma(i+1);
193 #ifdef R__HAS_MATHMORE
194 return ROOT::Math::noncentral_chisquared_pdf(_x,k,lambda);
196 coutF(Eval) <<
"RooNonCentralChisquare: ForceSum must be set" << endl;
204 Int_t RooNonCentralChiSquare::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars,
const char* )
const
206 if (matchArgs(allVars,analVars,x))
return 1 ;
212 Double_t RooNonCentralChiSquare::analyticalIntegral(Int_t code,
const char* rangeName)
const
216 Double_t xmin = x.min(rangeName);
217 Double_t xmax = x.max(rangeName);
223 return (ROOT::Math::chisquared_cdf(xmax,k) - ROOT::Math::chisquared_cdf(xmin,k));
234 double errorTol = fErrorTol;
235 int MaxIters = fMaxIters;
237 int iDominant = (int) TMath::Floor(lambda/2);
240 for(
int i = iDominant; ; ++i){
241 ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)
242 *( ROOT::Math::chisquared_cdf(xmax,k+2*i)/TMath::Gamma(i+1)
243 - ROOT::Math::chisquared_cdf(xmin,k+2*i)/TMath::Gamma(i+1) );
246 if(ithTerm/sum < errorTol)
249 if( i>iDominant+MaxIters){
250 if(!fHasIssuedConvWarning){
251 fHasIssuedConvWarning=
true;
252 coutW(Eval) <<
"RooNonCentralChiSquare Normalization did not converge: for k="<<k
253 <<
", lambda="<<lambda <<
" fractional error = " << ithTerm/sum
254 <<
"\n either adjust tolerance with SetErrorTolerance(tol) or max_iter with SetMaxIter(max_it)"
261 for(
int i = iDominant - 1; i >= 0; --i){
262 ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)
263 *( ROOT::Math::chisquared_cdf(xmax,k+2*i)/TMath::Gamma(i+1)
264 -ROOT::Math::chisquared_cdf(xmin,k+2*i)/TMath::Gamma(i+1) );