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 #include "Riostream.h"
00032
00033 #include "RooFit.h"
00034 #include "RooProjectedPdf.h"
00035 #include "RooMsgService.h"
00036 #include "RooAbsReal.h"
00037 #include "RooRealVar.h"
00038 #include "RooNameReg.h"
00039
00040
00041 ClassImp(RooProjectedPdf)
00042 ;
00043
00044
00045
00046 RooProjectedPdf::RooProjectedPdf() : _curNormSet(0)
00047 {
00048
00049 }
00050
00051
00052
00053
00054 RooProjectedPdf::RooProjectedPdf(const char *name, const char *title, RooAbsReal& _intpdf, const RooArgSet& intObs) :
00055 RooAbsPdf(name,title),
00056 intpdf("!IntegratedPdf","intpdf",this,_intpdf,kFALSE,kFALSE),
00057 intobs("!IntegrationObservables","intobs",this,kFALSE,kFALSE),
00058 deps("!Dependents","deps",this,kTRUE,kTRUE),
00059 _cacheMgr(this,10)
00060 {
00061
00062
00063 intobs.add(intObs) ;
00064
00065
00066 RooArgSet* tmpdeps = _intpdf.getParameters(intObs) ;
00067 deps.add(*tmpdeps) ;
00068 delete tmpdeps ;
00069 }
00070
00071
00072
00073
00074 RooProjectedPdf::RooProjectedPdf(const RooProjectedPdf& other, const char* name) :
00075 RooAbsPdf(other,name),
00076 intpdf("!IntegratedPdf",this,other.intpdf),
00077 intobs("!IntegrationObservable",this,other.intobs),
00078 deps("!Dependents",this,other.deps),
00079 _cacheMgr(other._cacheMgr,this)
00080 {
00081
00082 }
00083
00084
00085
00086
00087 Double_t RooProjectedPdf::getVal(const RooArgSet* set) const
00088 {
00089
00090
00091 _curNormSet = (RooArgSet*)set ;
00092
00093 return RooAbsPdf::getVal(set) ;
00094 }
00095
00096
00097
00098
00099 Double_t RooProjectedPdf::evaluate() const
00100 {
00101
00102
00103
00104 int code ;
00105 const RooAbsReal* proj = getProjection(&intobs,_curNormSet,0,code) ;
00106
00107 return proj->getVal() ;
00108 }
00109
00110
00111
00112
00113 const RooAbsReal* RooProjectedPdf::getProjection(const RooArgSet* iset, const RooArgSet* nset, const char* rangeName, int& code) const
00114 {
00115
00116
00117
00118
00119
00120
00121
00122 Int_t sterileIdx(-1) ;
00123 CacheElem* cache = (CacheElem*) _cacheMgr.getObj(iset,nset,&sterileIdx,RooNameReg::ptr(rangeName)) ;
00124 if (cache) {
00125 code = _cacheMgr.lastIndex() ;
00126 return static_cast<const RooAbsReal*>(cache->_projection);
00127 }
00128
00129 RooArgSet* nset2 = intpdf.arg().getObservables(*nset) ;
00130
00131 if (iset) {
00132 nset2->add(*iset) ;
00133 }
00134 RooAbsReal* proj = intpdf.arg().createIntegral(iset?*iset:RooArgSet(),nset2,0,rangeName) ;
00135 delete nset2 ;
00136
00137 cache = new CacheElem ;
00138 cache->_projection = proj ;
00139
00140 code = _cacheMgr.setObj(iset,nset,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ;
00141
00142 coutI(Integration) << "RooProjectedPdf::getProjection(" << GetName() << ") creating new projection " << proj->GetName() << " with code " << code << endl ;
00143
00144 return proj ;
00145 }
00146
00147
00148
00149
00150 RooAbsPdf* RooProjectedPdf::createProjection(const RooArgSet& iset)
00151 {
00152
00153
00154
00155
00156
00157 RooArgSet combiset(iset) ;
00158 combiset.add(intobs) ;
00159 return static_cast<RooAbsPdf&>( const_cast<RooAbsReal&>(intpdf.arg()) ).createProjection(combiset) ;
00160 }
00161
00162
00163
00164
00165 Bool_t RooProjectedPdf::forceAnalyticalInt(const RooAbsArg& ) const
00166 {
00167
00168
00169 return kTRUE ;
00170 }
00171
00172
00173
00174
00175 Int_t RooProjectedPdf::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName) const
00176 {
00177
00178
00179 analVars.add(allVars) ;
00180
00181
00182 int code ;
00183 RooArgSet allVars2(allVars) ;
00184 allVars2.add(intobs) ;
00185 getProjection(&allVars2,normSet,rangeName,code) ;
00186
00187 return code+1 ;
00188 }
00189
00190
00191
00192
00193 Double_t RooProjectedPdf::analyticalIntegralWN(Int_t code, const RooArgSet* , const char* rangeName) const
00194 {
00195
00196
00197 CacheElem *cache = (CacheElem*) _cacheMgr.getObjByIndex(code-1) ;
00198
00199 if (cache) {
00200 Double_t ret= cache->_projection->getVal() ;
00201 return ret ;
00202 } else {
00203
00204 RooArgSet* vars = getParameters(RooArgSet()) ;
00205 vars->add(intobs) ;
00206 RooArgSet* iset = _cacheMgr.nameSet1ByIndex(code-1)->select(*vars) ;
00207 RooArgSet* nset = _cacheMgr.nameSet2ByIndex(code-1)->select(*vars) ;
00208
00209 Int_t code2(-1) ;
00210 const RooAbsReal* proj = getProjection(iset,nset,rangeName,code2) ;
00211
00212 delete vars ;
00213 delete nset ;
00214 delete iset ;
00215
00216 Double_t ret = proj->getVal() ;
00217 return ret ;
00218 }
00219
00220 }
00221
00222
00223
00224
00225 Int_t RooProjectedPdf::getGenerator(const RooArgSet& , RooArgSet& , Bool_t ) const
00226 {
00227
00228 return 0 ;
00229 }
00230
00231
00232
00233
00234 void RooProjectedPdf::generateEvent(Int_t )
00235 {
00236
00237 return;
00238 }
00239
00240
00241
00242
00243 Bool_t RooProjectedPdf::redirectServersHook(const RooAbsCollection& newServerList, Bool_t ,
00244 Bool_t , Bool_t )
00245 {
00246
00247
00248
00249
00250
00251 RooAbsArg* newPdf = newServerList.find(intpdf.arg().GetName()) ;
00252 if (newPdf) {
00253
00254
00255 RooArgSet olddeps(deps) ;
00256 RooArgSet* newdeps = newPdf->getParameters(intobs) ;
00257 RooArgSet* common = (RooArgSet*) newdeps->selectCommon(deps) ;
00258 newdeps->remove(*common,kTRUE,kTRUE) ;
00259 olddeps.remove(*common,kTRUE,kTRUE) ;
00260
00261
00262 if (newdeps->getSize()>0) {
00263 deps.add(*newdeps) ;
00264 }
00265 if (olddeps.getSize()>0) {
00266 deps.remove(olddeps,kTRUE,kTRUE) ;
00267 }
00268
00269 delete common ;
00270 delete newdeps ;
00271 }
00272
00273 return kFALSE ;
00274 }
00275
00276
00277
00278
00279 RooArgList RooProjectedPdf::CacheElem::containedArgs(Action)
00280 {
00281
00282 RooArgList ret(*_projection) ;
00283 return ret ;
00284 }
00285
00286
00287
00288
00289 void RooProjectedPdf::printMetaArgs(ostream& os) const
00290 {
00291
00292
00293
00294 os << "Int " << intpdf.arg().GetName() ;
00295
00296 os << " d" ;
00297 os << intobs ;
00298 os << " " ;
00299
00300 }
00301
00302
00303
00304
00305
00306 void RooProjectedPdf::CacheElem::printCompactTreeHook(ostream& os, const char* indent, Int_t curElem, Int_t maxElem)
00307 {
00308
00309
00310 if (curElem==0) {
00311 os << indent << "RooProjectedPdf begin projection cache" << endl ;
00312 }
00313
00314 TString indent2(indent) ;
00315 indent2 += Form("[%d] ",curElem) ;
00316
00317 _projection->printCompactTree(os,indent2) ;
00318
00319 if(curElem==maxElem) {
00320 os << indent << "RooProjectedPdf end projection cache" << endl ;
00321 }
00322 }
00323
00324