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
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103 #include "Riostream.h"
00104
00105 #include "RooFit.h"
00106 #include "RooFFTConvPdf.h"
00107 #include "RooAbsReal.h"
00108 #include "RooMsgService.h"
00109 #include "RooDataHist.h"
00110 #include "RooHistPdf.h"
00111 #include "RooRealVar.h"
00112 #include "TComplex.h"
00113 #include "TVirtualFFT.h"
00114 #include "RooGenContext.h"
00115 #include "RooConvGenContext.h"
00116 #include "RooBinning.h"
00117 #include "RooLinearVar.h"
00118 #include "RooCustomizer.h"
00119 #include "RooGlobalFunc.h"
00120 #include "RooLinearVar.h"
00121 #include "RooConstVar.h"
00122 #include "TClass.h"
00123 #include "TSystem.h"
00124
00125 using namespace std ;
00126
00127 ClassImp(RooFFTConvPdf)
00128
00129
00130
00131
00132 RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooRealVar& convVar, RooAbsPdf& pdf1, RooAbsPdf& pdf2, Int_t ipOrder) :
00133 RooAbsCachedPdf(name,title,ipOrder),
00134 _x("!x","Convolution Variable",this,convVar),
00135 _xprime("!xprime","External Convolution Variable",this,0),
00136 _pdf1("!pdf1","pdf1",this,pdf1,kFALSE),
00137 _pdf2("!pdf2","pdf2",this,pdf2,kFALSE),
00138 _params("!params","effective parameters",this),
00139 _bufFrac(0.1),
00140 _bufStrat(Extend),
00141 _shift1(0),
00142 _shift2(0),
00143 _cacheObs("!cacheObs","Cached observables",this,kFALSE,kFALSE)
00144 {
00145
00146
00147
00148
00149 if (!convVar.hasBinning("cache")) {
00150 convVar.setBinning(convVar.getBinning(),"cache") ;
00151 }
00152
00153 _shift2 = (convVar.getMax()+convVar.getMin())/2 ;
00154
00155 calcParams() ;
00156
00157 }
00158
00159
00160
00161
00162 RooFFTConvPdf::RooFFTConvPdf(const char *name, const char *title, RooAbsReal& pdfConvVar, RooRealVar& convVar, RooAbsPdf& pdf1, RooAbsPdf& pdf2, Int_t ipOrder) :
00163 RooAbsCachedPdf(name,title,ipOrder),
00164 _x("!x","Convolution Variable",this,convVar,kFALSE,kFALSE),
00165 _xprime("!xprime","External Convolution Variable",this,pdfConvVar),
00166 _pdf1("!pdf1","pdf1",this,pdf1,kFALSE),
00167 _pdf2("!pdf2","pdf2",this,pdf2,kFALSE),
00168 _params("!params","effective parameters",this),
00169 _bufFrac(0.1),
00170 _bufStrat(Extend),
00171 _shift1(0),
00172 _shift2(0),
00173 _cacheObs("!cacheObs","Cached observables",this,kFALSE,kFALSE)
00174 {
00175
00176
00177
00178
00179 if (!convVar.hasBinning("cache")) {
00180 convVar.setBinning(convVar.getBinning(),"cache") ;
00181 }
00182
00183 _shift2 = (convVar.getMax()+convVar.getMin())/2 ;
00184
00185 calcParams() ;
00186 }
00187
00188
00189
00190
00191 RooFFTConvPdf::RooFFTConvPdf(const RooFFTConvPdf& other, const char* name) :
00192 RooAbsCachedPdf(other,name),
00193 _x("!x",this,other._x),
00194 _xprime("!xprime",this,other._xprime),
00195 _pdf1("!pdf1",this,other._pdf1),
00196 _pdf2("!pdf2",this,other._pdf2),
00197 _params("!params",this,other._params),
00198 _bufFrac(other._bufFrac),
00199 _bufStrat(other._bufStrat),
00200 _shift1(other._shift1),
00201 _shift2(other._shift2),
00202 _cacheObs("!cacheObs",this,other._cacheObs)
00203 {
00204
00205 }
00206
00207
00208
00209
00210 RooFFTConvPdf::~RooFFTConvPdf()
00211 {
00212
00213 }
00214
00215
00216
00217
00218 const char* RooFFTConvPdf::inputBaseName() const
00219 {
00220
00221
00222 static TString name ;
00223 name = _pdf1.arg().GetName() ;
00224 name.Append("_CONV_") ;
00225 name.Append(_pdf2.arg().GetName()) ;
00226 return name.Data() ;
00227 }
00228
00229
00230
00231
00232
00233 RooFFTConvPdf::PdfCacheElem* RooFFTConvPdf::createCache(const RooArgSet* nset) const
00234 {
00235
00236 return new FFTCacheElem(*this,nset) ;
00237 }
00238
00239
00240
00241
00242
00243 RooFFTConvPdf::FFTCacheElem::FFTCacheElem(const RooFFTConvPdf& self, const RooArgSet* nsetIn) :
00244 PdfCacheElem(self,nsetIn),
00245 fftr2c1(0),fftr2c2(0),fftc2r(0)
00246 {
00247
00248
00249 RooAbsPdf* clonePdf1 = (RooAbsPdf*) self._pdf1.arg().cloneTree() ;
00250 RooAbsPdf* clonePdf2 = (RooAbsPdf*) self._pdf2.arg().cloneTree() ;
00251 clonePdf1->attachDataSet(*hist()) ;
00252 clonePdf2->attachDataSet(*hist()) ;
00253
00254
00255 RooRealVar* convObs = (RooRealVar*) hist()->get()->find(self._x.arg().GetName()) ;
00256
00257
00258 string refName = Form("refrange_fft_%s",self.GetName()) ;
00259 convObs->setRange(refName.c_str(),convObs->getMin(),convObs->getMax()) ;
00260
00261 if (self._shift1!=0) {
00262 RooLinearVar* shiftObs1 = new RooLinearVar(Form("%s_shifted_FFTBuffer1",convObs->GetName()),"shiftObs1",
00263 *convObs,RooFit::RooConst(1),RooFit::RooConst(-1*self._shift1)) ;
00264
00265 RooArgSet clonedBranches1 ;
00266 RooCustomizer cust(*clonePdf1,"fft") ;
00267 cust.replaceArg(*convObs,*shiftObs1) ;
00268
00269 pdf1Clone = (RooAbsPdf*) cust.build() ;
00270
00271 pdf1Clone->addOwnedComponents(*shiftObs1) ;
00272 pdf1Clone->addOwnedComponents(*clonePdf1) ;
00273
00274 } else {
00275 pdf1Clone = clonePdf1 ;
00276 }
00277
00278 if (self._shift2!=0) {
00279 RooLinearVar* shiftObs2 = new RooLinearVar(Form("%s_shifted_FFTBuffer2",convObs->GetName()),"shiftObs2",
00280 *convObs,RooFit::RooConst(1),RooFit::RooConst(-1*self._shift2)) ;
00281
00282 RooArgSet clonedBranches2 ;
00283 RooCustomizer cust(*clonePdf2,"fft") ;
00284 cust.replaceArg(*convObs,*shiftObs2) ;
00285
00286 pdf1Clone->addOwnedComponents(*shiftObs2) ;
00287 pdf1Clone->addOwnedComponents(*clonePdf2) ;
00288
00289 pdf2Clone = (RooAbsPdf*) cust.build() ;
00290
00291 } else {
00292 pdf2Clone = clonePdf2 ;
00293 }
00294
00295
00296
00297 RooArgSet* fftParams = self.getParameters(*convObs) ;
00298
00299
00300
00301 fftParams->remove(*hist()->get(),kTRUE,kTRUE) ;
00302
00303 pdf1Clone->recursiveRedirectServers(*fftParams) ;
00304 pdf2Clone->recursiveRedirectServers(*fftParams) ;
00305 pdf1Clone->fixAddCoefRange(refName.c_str()) ;
00306 pdf2Clone->fixAddCoefRange(refName.c_str()) ;
00307
00308 delete fftParams ;
00309
00310
00311
00312
00313 Int_t N = convObs->numBins() ;
00314 Int_t Nbuf = static_cast<Int_t>((N*self.bufferFraction())/2 + 0.5) ;
00315 Double_t obw = (convObs->getMax() - convObs->getMin())/N ;
00316 Int_t N2 = N+2*Nbuf ;
00317
00318 scanBinning = new RooUniformBinning (convObs->getMin()-Nbuf*obw,convObs->getMax()+Nbuf*obw,N2) ;
00319 histBinning = convObs->getBinning().clone() ;
00320
00321
00322
00323 hist()->setDirtyProp(kFALSE) ;
00324 convObs->setOperMode(ADirty,kTRUE) ;
00325 }
00326
00327
00328
00329 TString RooFFTConvPdf::histNameSuffix() const
00330 {
00331
00332 return TString(Form("_BufFrac%3.1f_BufStrat%d",_bufFrac,_bufStrat)) ;
00333 }
00334
00335
00336
00337
00338 RooArgList RooFFTConvPdf::FFTCacheElem::containedArgs(Action a)
00339 {
00340
00341 RooArgList ret(PdfCacheElem::containedArgs(a)) ;
00342
00343 ret.add(*pdf1Clone) ;
00344 ret.add(*pdf2Clone) ;
00345 if (pdf1Clone->ownedComponents()) {
00346 ret.add(*pdf1Clone->ownedComponents()) ;
00347 }
00348 if (pdf2Clone->ownedComponents()) {
00349 ret.add(*pdf2Clone->ownedComponents()) ;
00350 }
00351
00352 return ret ;
00353 }
00354
00355
00356
00357 RooFFTConvPdf::FFTCacheElem::~FFTCacheElem()
00358 {
00359 delete fftr2c1 ;
00360 delete fftr2c2 ;
00361 delete fftc2r ;
00362
00363 delete pdf1Clone ;
00364 delete pdf2Clone ;
00365
00366 delete histBinning ;
00367 delete scanBinning ;
00368
00369 }
00370
00371
00372
00373
00374
00375 void RooFFTConvPdf::fillCacheObject(RooAbsCachedPdf::PdfCacheElem& cache) const
00376 {
00377
00378 RooDataHist& cacheHist = *cache.hist() ;
00379
00380 ((FFTCacheElem&)cache).pdf1Clone->setOperMode(ADirty,kTRUE) ;
00381 ((FFTCacheElem&)cache).pdf2Clone->setOperMode(ADirty,kTRUE) ;
00382
00383
00384 RooArgSet otherObs ;
00385 RooArgSet(*cacheHist.get()).snapshot(otherObs) ;
00386
00387 RooAbsArg* histArg = otherObs.find(_x.arg().GetName()) ;
00388 if (histArg) {
00389 otherObs.remove(*histArg,kTRUE,kTRUE) ;
00390 delete histArg ;
00391 }
00392
00393
00394 if (otherObs.getSize()==0) {
00395 fillCacheSlice((FFTCacheElem&)cache,RooArgSet()) ;
00396 return ;
00397 }
00398
00399
00400
00401
00402
00403 Int_t n = otherObs.getSize() ;
00404 Int_t* binCur = new Int_t[n+1] ;
00405 Int_t* binMax = new Int_t[n+1] ;
00406 Int_t curObs = 0 ;
00407
00408 RooAbsLValue** obsLV = new RooAbsLValue*[n] ;
00409 TIterator* iter = otherObs.createIterator() ;
00410 RooAbsArg* arg ;
00411 Int_t i(0) ;
00412 while((arg=(RooAbsArg*)iter->Next())) {
00413 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(arg) ;
00414 obsLV[i] = lvarg ;
00415 binCur[i] = 0 ;
00416
00417 binMax[i] = lvarg->numBins(binningName())-1 ;
00418 i++ ;
00419 }
00420 delete iter ;
00421
00422 Bool_t loop(kTRUE) ;
00423 while(loop) {
00424
00425 for (Int_t j=0 ; j<n ; j++) { obsLV[j]->setBin(binCur[j],binningName()) ; }
00426
00427
00428
00429
00430 fillCacheSlice((FFTCacheElem&)cache,otherObs) ;
00431
00432
00433 while(binCur[curObs]==binMax[curObs]) {
00434
00435
00436 binCur[curObs]=0 ;
00437 curObs++ ;
00438
00439
00440 if (curObs==n) {
00441 loop=kFALSE ;
00442 break ;
00443 }
00444 }
00445
00446
00447 binCur[curObs]++ ;
00448 curObs=0 ;
00449
00450 }
00451
00452 delete[] obsLV ;
00453 delete[] binMax ;
00454 delete[] binCur ;
00455
00456 }
00457
00458
00459
00460 void RooFFTConvPdf::fillCacheSlice(FFTCacheElem& aux, const RooArgSet& slicePos) const
00461 {
00462
00463
00464
00465 RooDataHist& cacheHist = *aux.hist() ;
00466
00467
00468
00469
00470
00471
00472
00473
00474
00475
00476
00477
00478 Int_t N,N2,binShift1,binShift2 ;
00479
00480 RooRealVar* histX = (RooRealVar*) cacheHist.get()->find(_x.arg().GetName()) ;
00481 if (_bufStrat==Extend) histX->setBinning(*aux.scanBinning) ;
00482 Double_t* input1 = scanPdf((RooRealVar&)_x.arg(),*aux.pdf1Clone,cacheHist,slicePos,N,N2,binShift1,_shift1) ;
00483 Double_t* input2 = scanPdf((RooRealVar&)_x.arg(),*aux.pdf2Clone,cacheHist,slicePos,N,N2,binShift2,_shift2) ;
00484 if (_bufStrat==Extend) histX->setBinning(*aux.histBinning) ;
00485
00486
00487 if (!aux.fftr2c1) {
00488 aux.fftr2c1 = TVirtualFFT::FFT(1, &N2, "R2CK");
00489 aux.fftr2c2 = TVirtualFFT::FFT(1, &N2, "R2CK");
00490 aux.fftc2r = TVirtualFFT::FFT(1, &N2, "C2RK");
00491 }
00492
00493
00494 aux.fftr2c1->SetPoints(input1);
00495 aux.fftr2c1->Transform();
00496
00497
00498 aux.fftr2c2->SetPoints(input2);
00499 aux.fftr2c2->Transform();
00500
00501
00502
00503 for (Int_t i=0 ; i<N2/2+1 ; i++) {
00504 Double_t re1,re2,im1,im2 ;
00505 aux.fftr2c1->GetPointComplex(i,re1,im1) ;
00506 aux.fftr2c2->GetPointComplex(i,re2,im2) ;
00507 Double_t re = re1*re2 - im1*im2 ;
00508 Double_t im = re1*im2 + re2*im1 ;
00509 TComplex t(re,im) ;
00510 aux.fftc2r->SetPointComplex(i,t) ;
00511 }
00512
00513
00514 aux.fftc2r->Transform() ;
00515
00516 Int_t totalShift = binShift1 + (N2-N)/2 ;
00517
00518
00519 TIterator* iter = const_cast<RooDataHist&>(cacheHist).sliceIterator(const_cast<RooAbsReal&>(_x.arg()),slicePos) ;
00520 for (Int_t i =0 ; i<N ; i++) {
00521
00522
00523 Int_t j = i + totalShift ;
00524 while (j<0) j+= N2 ;
00525 while (j>=N2) j-= N2 ;
00526
00527 iter->Next() ;
00528 cacheHist.set(aux.fftc2r->GetPointReal(j)) ;
00529 }
00530 delete iter ;
00531
00532
00533
00534
00535 delete[] input1 ;
00536 delete[] input2 ;
00537
00538 }
00539
00540
00541
00542 Double_t* RooFFTConvPdf::scanPdf(RooRealVar& obs, RooAbsPdf& pdf, const RooDataHist& hist, const RooArgSet& slicePos,
00543 Int_t& N, Int_t& N2, Int_t& zeroBin, Double_t shift) const
00544 {
00545
00546
00547
00548
00549
00550
00551 RooRealVar* histX = (RooRealVar*) hist.get()->find(obs.GetName()) ;
00552
00553
00554 N = histX->numBins(binningName()) ;
00555 Int_t Nbuf = static_cast<Int_t>((N*bufferFraction())/2 + 0.5) ;
00556 N2 = N+2*Nbuf ;
00557
00558
00559
00560 Double_t* array = new Double_t[N2] ;
00561
00562
00563 hist.get(slicePos) ;
00564
00565
00566 zeroBin = 0 ;
00567 if (histX->getMax()>=0 && histX->getMin()<=0) {
00568 zeroBin = histX->getBinning().binNumber(0) ;
00569 } else if (histX->getMin()>0) {
00570 Double_t bw = (histX->getMax() - histX->getMin())/N2 ;
00571 zeroBin = Int_t(-histX->getMin()/bw) ;
00572 } else {
00573 Double_t bw = (histX->getMax() - histX->getMin())/N2 ;
00574 zeroBin = Int_t(-1*histX->getMax()/bw) ;
00575 }
00576
00577 Int_t binShift = Int_t((N2* shift) / (histX->getMax()-histX->getMin())) ;
00578
00579 zeroBin += binShift ;
00580 while(zeroBin>=N2) zeroBin-= N2 ;
00581 while(zeroBin<0) zeroBin+= N2 ;
00582
00583
00584 Double_t *tmp = new Double_t[N2] ;
00585 Int_t k(0) ;
00586 switch(_bufStrat) {
00587
00588 case Extend:
00589
00590 for (k=0 ; k<N2 ; k++) {
00591 histX->setBin(k) ;
00592 tmp[k] = pdf.getVal(hist.get()) ;
00593 }
00594 break ;
00595
00596 case Flat:
00597
00598
00599 {
00600 histX->setBin(0) ;
00601 Double_t val = pdf.getVal(hist.get()) ;
00602 for (k=0 ; k<Nbuf ; k++) {
00603 tmp[k] = val ;
00604 }
00605 for (k=0 ; k<N ; k++) {
00606 histX->setBin(k) ;
00607 tmp[k+Nbuf] = pdf.getVal(hist.get()) ;
00608 }
00609 histX->setBin(N-1) ;
00610 val = pdf.getVal(hist.get()) ;
00611 for (k=0 ; k<Nbuf ; k++) {
00612 tmp[N+Nbuf+k] = val ;
00613 }
00614 }
00615 break ;
00616
00617 case Mirror:
00618
00619
00620 for (k=0 ; k<N ; k++) {
00621 histX->setBin(k) ;
00622 tmp[k+Nbuf] = pdf.getVal(hist.get()) ;
00623 }
00624 for (k=1 ; k<=Nbuf ; k++) {
00625 histX->setBin(k) ;
00626 tmp[Nbuf-k] = pdf.getVal(hist.get()) ;
00627 histX->setBin(N-k) ;
00628 tmp[Nbuf+N+k-1] = pdf.getVal(hist.get()) ;
00629 }
00630 break ;
00631 }
00632
00633
00634 for (Int_t i=0 ; i<N2 ; i++) {
00635
00636 Int_t j = i - (zeroBin) ;
00637 if (j<0) j+= N2 ;
00638 if (j>=N2) j-= N2 ;
00639 array[i] = tmp[j] ;
00640 }
00641
00642
00643 delete[] tmp ;
00644 return array ;
00645 }
00646
00647
00648
00649
00650 RooArgSet* RooFFTConvPdf::actualObservables(const RooArgSet& nset) const
00651 {
00652
00653
00654
00655
00656
00657
00658
00659
00660
00661
00662
00663
00664
00665 RooArgSet* obs1 = _pdf1.arg().getObservables(nset) ;
00666 RooArgSet* obs2 = _pdf2.arg().getObservables(nset) ;
00667 obs1->add(*obs2,kTRUE) ;
00668
00669
00670 if (nset.contains(_x.arg())) {
00671
00672
00673 TIterator* iter = obs1->createIterator() ;
00674 RooAbsArg* arg ;
00675 RooArgSet killList ;
00676 while((arg=(RooAbsArg*)iter->Next())) {
00677 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
00678 killList.add(*arg) ;
00679 }
00680 }
00681 delete iter ;
00682 obs1->remove(killList) ;
00683
00684
00685 obs1->add(_x.arg(),kTRUE) ;
00686 delete obs2 ;
00687
00688 } else {
00689
00690
00691 if (_cacheObs.getSize()>0) {
00692 TIterator* iter = obs1->createIterator() ;
00693 RooAbsArg* arg ;
00694 RooArgSet killList ;
00695 while((arg=(RooAbsArg*)iter->Next())) {
00696 if (arg->IsA()->InheritsFrom(RooAbsReal::Class()) && !_cacheObs.find(arg->GetName())) {
00697 killList.add(*arg) ;
00698 }
00699 }
00700 delete iter ;
00701 obs1->remove(killList) ;
00702 }
00703
00704
00705
00706 obs1->add(_x.arg(),kTRUE) ;
00707 delete obs2 ;
00708
00709 }
00710
00711 return obs1 ;
00712 }
00713
00714
00715
00716
00717 RooArgSet* RooFFTConvPdf::actualParameters(const RooArgSet& nset) const
00718 {
00719
00720
00721
00722
00723 RooArgSet* vars = getVariables() ;
00724 RooArgSet* obs = actualObservables(nset) ;
00725 vars->remove(*obs) ;
00726 delete obs ;
00727
00728 return vars ;
00729 }
00730
00731
00732
00733
00734 RooAbsArg& RooFFTConvPdf::pdfObservable(RooAbsArg& histObservable) const
00735 {
00736
00737
00738
00739 if (_xprime.absArg() && string(histObservable.GetName())==_x.absArg()->GetName()) {
00740 return (*_xprime.absArg()) ;
00741 }
00742 return histObservable ;
00743 }
00744
00745
00746
00747
00748 RooAbsGenContext* RooFFTConvPdf::genContext(const RooArgSet &vars, const RooDataSet *prototype,
00749 const RooArgSet* auxProto, Bool_t verbose) const
00750 {
00751
00752
00753
00754
00755
00756
00757 RooArgSet vars2(vars) ;
00758 vars2.remove(_x.arg(),kTRUE,kTRUE) ;
00759 Int_t numAddDep = vars2.getSize() ;
00760
00761 RooArgSet dummy ;
00762 Bool_t pdfCanDir = (((RooAbsPdf&)_pdf1.arg()).getGenerator(_x.arg(),dummy) != 0 && \
00763 ((RooAbsPdf&)_pdf1.arg()).isDirectGenSafe(_x.arg())) ;
00764 Bool_t resCanDir = (((RooAbsPdf&)_pdf2.arg()).getGenerator(_x.arg(),dummy) !=0 &&
00765 ((RooAbsPdf&)_pdf2.arg()).isDirectGenSafe(_x.arg())) ;
00766
00767 if (pdfCanDir) {
00768 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f " << _pdf1.arg().GetName()
00769 << " has internal generator that is safe to use in current context" << endl ;
00770 }
00771 if (resCanDir) {
00772 cxcoutI(Generation) << "RooFFTConvPdf::genContext() input p.d.f. " << _pdf2.arg().GetName()
00773 << " has internal generator that is safe to use in current context" << endl ;
00774 }
00775 if (numAddDep>0) {
00776 cxcoutI(Generation) << "RooFFTConvPdf::genContext() generation requested for observables other than the convolution observable " << _x.arg().GetName() << endl ;
00777 }
00778
00779
00780 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
00781
00782
00783 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting accept/reject generator context because one or both of the input "
00784 << "p.d.f.s cannot use internal generator and/or "
00785 << "observables other than the convolution variable are requested for generation" << endl ;
00786 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
00787 }
00788
00789
00790 cxcoutI(Generation) << "RooFFTConvPdf::genContext() selecting specialized convolution generator context as both input "
00791 << "p.d.fs are safe for internal generator and only "
00792 << "the convolution observables is requested for generation" << endl ;
00793 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
00794 }
00795
00796
00797
00798
00799 void RooFFTConvPdf::setBufferFraction(Double_t frac)
00800 {
00801
00802
00803
00804 if (frac<0) {
00805 coutE(InputArguments) << "RooFFTConvPdf::setBufferFraction(" << GetName() << ") fraction should be greater than or equal to zero" << endl ;
00806 return ;
00807 }
00808 _bufFrac = frac ;
00809
00810
00811 _cacheMgr.sterilize() ;
00812 }
00813
00814
00815
00816 void RooFFTConvPdf::setBufferStrategy(BufStrat bs)
00817 {
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827
00828 _bufStrat = bs ;
00829 }
00830
00831
00832
00833
00834 void RooFFTConvPdf::printMetaArgs(ostream& os) const
00835 {
00836
00837
00838
00839 os << _pdf1.arg().GetName() << "(" << _x.arg().GetName() << ") (*) " << _pdf2.arg().GetName() << "(" << _x.arg().GetName() << ") " ;
00840 }
00841
00842
00843
00844
00845 void RooFFTConvPdf::calcParams()
00846 {
00847
00848
00849 RooArgSet* params1 = _pdf1.arg().getParameters(_x.arg()) ;
00850 RooArgSet* params2 = _pdf2.arg().getParameters(_x.arg()) ;
00851 _params.removeAll() ;
00852 _params.add(*params1) ;
00853 _params.add(*params2,kTRUE) ;
00854 delete params1 ;
00855 delete params2 ;
00856 }
00857
00858
00859
00860
00861 Bool_t RooFFTConvPdf::redirectServersHook(const RooAbsCollection& , Bool_t , Bool_t , Bool_t )
00862 {
00863
00864 return kFALSE ;
00865 }