00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044 #include "RooFit.h"
00045
00046 #include "Riostream.h"
00047 #include "TArrayD.h"
00048 #include <math.h>
00049
00050 #include "RooParametricStepFunction.h"
00051 #include "RooAbsReal.h"
00052 #include "RooRealVar.h"
00053 #include "RooArgList.h"
00054
00055 ClassImp(RooParametricStepFunction)
00056 ;
00057
00058
00059
00060 RooParametricStepFunction::RooParametricStepFunction(const char* name, const char* title,
00061 RooAbsReal& x, const RooArgList& coefList, TArrayD& limits, Int_t nBins) :
00062 RooAbsPdf(name, title),
00063 _x("x", "Dependent", this, x),
00064 _coefList("coefList","List of coefficients",this),
00065 _nBins(nBins)
00066 {
00067
00068 _coefIter = _coefList.createIterator() ;
00069
00070
00071 if (_nBins<0) {
00072 cout << "RooParametricStepFunction::ctor(" << GetName()
00073 << ") WARNING: nBins must be >=0, setting value to 0" << endl ;
00074 _nBins=0 ;
00075 }
00076
00077 TIterator* coefIter = coefList.createIterator() ;
00078 RooAbsArg* coef ;
00079 while((coef = (RooAbsArg*)coefIter->Next())) {
00080 if (!dynamic_cast<RooAbsReal*>(coef)) {
00081 cout << "RooParametricStepFunction::ctor(" << GetName() << ") ERROR: coefficient " << coef->GetName()
00082 << " is not of type RooAbsReal" << endl ;
00083 assert(0) ;
00084 }
00085 _coefList.add(*coef) ;
00086 }
00087 delete coefIter ;
00088
00089
00090 limits.Copy(_limits);
00091
00092 }
00093
00094
00095
00096
00097 RooParametricStepFunction::RooParametricStepFunction(const RooParametricStepFunction& other, const char* name) :
00098 RooAbsPdf(other, name),
00099 _x("x", this, other._x),
00100 _coefList("coefList",this,other._coefList),
00101 _nBins(other._nBins)
00102 {
00103
00104 _coefIter = _coefList.createIterator();
00105 (other._limits).Copy(_limits);
00106 }
00107
00108
00109
00110
00111 RooParametricStepFunction::~RooParametricStepFunction()
00112 {
00113
00114 delete _coefIter ;
00115 }
00116
00117
00118
00119
00120 Int_t RooParametricStepFunction::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* ) const
00121 {
00122 if (matchArgs(allVars, analVars, _x)) return 1;
00123 return 0;
00124 }
00125
00126
00127
00128
00129 Double_t RooParametricStepFunction::analyticalIntegral(Int_t code, const char* rangeName) const
00130 {
00131 assert(code==1) ;
00132
00133
00134 if (!rangeName) {
00135 return 1.0 ;
00136 }
00137
00138
00139 Double_t xmin = _x.min(rangeName) ;
00140 Double_t xmax = _x.max(rangeName) ;
00141
00142 Double_t sum=0 ;
00143 Int_t i ;
00144 for (i=1 ; i<=_nBins ; i++) {
00145 Double_t binVal = (i<_nBins) ? (static_cast<RooRealVar*>(_coefList.at(i-1))->getVal()) : lastBinValue() ;
00146 if (_limits[i-1]>=xmin && _limits[i]<=xmax) {
00147
00148 sum += (_limits[i]-_limits[i-1])*binVal ;
00149 } else if (_limits[i-1]<xmin && _limits[i]>xmax) {
00150
00151 sum += (xmax-xmin)*binVal ;
00152
00153 return sum ;
00154 } else if (_limits[i-1]<xmin && _limits[i]<=xmax && _limits[i]>xmin) {
00155
00156 sum += (_limits[i]-xmin)*binVal ;
00157 } else if (_limits[i-1]>=xmin && _limits[i]>xmax && _limits[i-1]<xmax) {
00158 sum += (xmax-_limits[i-1])*binVal ;
00159
00160
00161 return sum ;
00162 }
00163 }
00164
00165 return sum;
00166 }
00167
00168
00169
00170
00171 Double_t RooParametricStepFunction::lastBinValue() const
00172 {
00173 Double_t sum(0.);
00174 Double_t binSize(0.);
00175 for (Int_t j=1;j<_nBins;j++){
00176 RooRealVar* tmp = (RooRealVar*) _coefList.at(j-1);
00177 binSize = _limits[j] - _limits[j-1];
00178 sum = sum + tmp->getVal()*binSize;
00179 }
00180 binSize = _limits[_nBins] - _limits[_nBins-1];
00181 return (1.0 - sum)/binSize;
00182 }
00183
00184
00185
00186
00187 Double_t RooParametricStepFunction::evaluate() const
00188 {
00189 Double_t xval(0.);
00190 xval = _x;
00191 Double_t value(0.);
00192 if (_x >= _limits[0] && _x < _limits[_nBins]){
00193
00194 for (Int_t i=1;i<=_nBins;i++){
00195 if (_x < _limits[i]){
00196
00197 if (i<_nBins) {
00198
00199 RooRealVar* tmp = (RooRealVar*) _coefList.at(i-1);
00200 value = tmp->getVal();
00201 break;
00202 } else {
00203
00204 Double_t sum(0.);
00205 Double_t binSize(0.);
00206 for (Int_t j=1;j<_nBins;j++){
00207 RooRealVar* tmp = (RooRealVar*) _coefList.at(j-1);
00208 binSize = _limits[j] - _limits[j-1];
00209 sum = sum + tmp->getVal()*binSize;
00210 }
00211 binSize = _limits[_nBins] - _limits[_nBins-1];
00212 value = (1.0 - sum)/binSize;
00213 if (value<=0.0){
00214 value = 0.000001;
00215
00216 }
00217 break;
00218 }
00219 }
00220 }
00221
00222 }
00223 return value;
00224
00225 }
00226
00227
00228
00229 Int_t RooParametricStepFunction::getnBins(){
00230 return _nBins;
00231 }
00232
00233
00234
00235 Double_t* RooParametricStepFunction::getLimits(){
00236 Double_t* limoutput = _limits.GetArray();
00237 return limoutput;
00238 }